LCOV - code coverage report
Current view: top level - ugbase/lib_disc/operator/linear_operator/multi_grid_solver - mg_solver.h (source / functions) Coverage Total Hit
Test: coverage.info Lines: 0.0 % 56 0
Test Date: 2025-09-21 23:31:46 Functions: 0.0 % 252 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__MULTI_GRID_SOLVER__MG_SOLVER__
      34              : #define __H__UG__LIB_DISC__MULTI_GRID_SOLVER__MG_SOLVER__
      35              : 
      36              : // extern includes
      37              : #include <vector>
      38              : #include <iostream>
      39              : 
      40              : // other ug4 modules
      41              : #include "common/common.h"
      42              : #ifdef UG_PARALLEL
      43              :         #include "lib_algebra/parallelization/parallel_storage_type.h"
      44              :         #include "lib_grid/parallelization/distributed_grid.h"
      45              : #endif
      46              : #include "lib_algebra/operator/interface/linear_iterator.h"
      47              : #include "lib_algebra/operator/interface/operator_inverse.h"
      48              : #include "lib_algebra/operator/interface/operator.h"
      49              : #include "lib_algebra/operator/preconditioner/jacobi.h"
      50              : #include "lib_algebra/operator/linear_solver/lu.h"
      51              : #include "lib_disc/dof_manager/dof_distribution.h"
      52              : #include "lib_disc/operator/linear_operator/transfer_interface.h"
      53              : //only for debugging!!!
      54              : #include "lib_grid/algorithms/debug_util.h"
      55              : 
      56              : // library intern headers
      57              : #include "lib_disc/function_spaces/grid_function_util.h"
      58              : #include "lib_disc/operator/linear_operator/assembled_linear_operator.h"
      59              : 
      60              : #include "mg_stats.h"
      61              : 
      62              : namespace ug{
      63              : 
      64              : ////////////////////////////////////////////////////////////////////////////////
      65              : // AssembledMultiGridCycle
      66              : ////////////////////////////////////////////////////////////////////////////////
      67              : 
      68              : /// geometric multi grid preconditioner
      69              : /**
      70              :  * This class implements one step of the geometric multi grid as a
      71              :  * preconditioner for linear iteration schemes such as linear iteration, cg
      72              :  * or bicgstab.
      73              :  *
      74              :  * The coarse grid spaces are build up according to the Approximation Space
      75              :  * that is set from outside. In addition an Assembling routine must be
      76              :  * specified that is used to assemble the coarse grid matrices.
      77              :  *
      78              :  * \tparam              TApproximationSpace             Type of Approximation Space
      79              :  * \tparam              TAlgebra                                Type of Algebra
      80              :  */
      81              : template <typename TDomain, typename TAlgebra>
      82              : class AssembledMultiGridCycle :
      83              :  public ILinearIterator<     typename TAlgebra::vector_type>
      84              : {
      85              :         public:
      86              :         ///     Domain
      87              :                 typedef TDomain domain_type;
      88              : 
      89              :         ///     Algebra type
      90              :                 typedef TAlgebra algebra_type;
      91              : 
      92              :         ///     Vector type
      93              :                 typedef typename algebra_type::vector_type vector_type;
      94              : 
      95              :         ///     Matrix type
      96              :                 typedef typename algebra_type::matrix_type matrix_type;
      97              : 
      98              :         ///     Grid Function type
      99              :                 typedef GridFunction<TDomain, TAlgebra> GF;
     100              : 
     101              :         ///     MGStats type
     102              :                 typedef MGStats<TDomain, TAlgebra> mg_stats_type;
     103              : 
     104              :         ///////////////////////////////////////////////////////////////////////////
     105              :         //      Setup
     106              :         ///////////////////////////////////////////////////////////////////////////
     107              : 
     108              :         public:
     109              :         /// constructor without arguments
     110              :                 AssembledMultiGridCycle();
     111              : 
     112              :         /// constructor setting approximation space
     113              :                 AssembledMultiGridCycle(SmartPtr<ApproximationSpace<TDomain> > approxSpace);
     114              : 
     115              :         ///     sets MGStats, an object which records statistics for individual iterations
     116              :         /**     Setting MGStats is optional. If set, a runtime overhead is incurred. Only
     117              :          * useful for debugging.
     118              :          * \sa MGStats*/
     119            0 :                 void set_mg_stats(SmartPtr<mg_stats_type> mgstats)
     120            0 :                         {m_mgstats = mgstats;}
     121              : 
     122              :         /// sets the approximation space
     123              :                 void set_approximation_space(SmartPtr<ApproximationSpace<TDomain> > approxSpace);
     124              : 
     125              :         /// sets the assembling procedure that is used to compute coarse grid matrices
     126            0 :                 void set_discretization(SmartPtr<IAssemble<TAlgebra> > spAss)
     127            0 :                         {m_spAss = spAss;}
     128              : 
     129              :         ///     sets the level where exact solving is performed in the mg cycle
     130            0 :                 void set_base_level(int baseLevel) {m_baseLev = baseLevel;}
     131              : 
     132              :         ///     sets the surface level (default is top-surface)
     133            0 :                 void set_surface_level(int topLevel) {m_surfaceLev = topLevel;}
     134              : 
     135              :         ///     sets the base solver that is used
     136            0 :                 void set_base_solver(SmartPtr<ILinearOperatorInverse<vector_type> > baseSolver)
     137            0 :                         {m_spBaseSolver = baseSolver;}
     138              : 
     139              :         ///     sets if the base solver is applied in parallel
     140            0 :                 void set_gathered_base_solver_if_ambiguous(bool bGathered) {m_bGatheredBaseIfAmbiguous = bGathered;}
     141              : 
     142              :         ///     sets if copies should be used to emulate a full-refined grid
     143            0 :                 void set_emulate_full_refined_grid(bool bEmulate){
     144            0 :                         if(bEmulate) m_GridLevelType = GridLevel::SURFACE;
     145            0 :                         else m_GridLevelType = GridLevel::LEVEL;
     146            0 :                 }
     147              : 
     148              :         ///     sets if RAP - Product used to build coarse grid matrices
     149            0 :                 void set_rap(bool bRAP) {m_bUseRAP = bRAP;}
     150              : 
     151              :         ///     sets if smoothing is performed on surface rim
     152            0 :                 void set_smooth_on_surface_rim(bool bSmooth) {m_bSmoothOnSurfaceRim = bSmooth;}
     153              : 
     154              :         ///     sets the cycle type (1 = V-cycle, 2 = W-cycle, ...)
     155            0 :                 void set_cycle_type(int type) {m_cycleType = type;}
     156              : 
     157              :         ///     sets the cycle type (1 = V-cycle, 2 = W-cycle, ...)
     158            0 :                 void set_cycle_type(const std::string& type) {
     159            0 :                         if(TrimString(type) == "V") {m_cycleType = _V_;}
     160            0 :                         else if(TrimString(type) == "W") {m_cycleType = _W_;}
     161            0 :                         else if(TrimString(type) == "F") {m_cycleType = _F_;}
     162            0 :                         else {UG_THROW("GMG::set_cycle_type: option '"<<type<<"' not supported.");}
     163            0 :                 }
     164              : 
     165              :         ///     sets if communication and computation should be overlaped
     166            0 :                 void set_comm_comp_overlap(bool bOverlap) {m_bCommCompOverlap = bOverlap;}
     167              : 
     168              :         ///     sets the number of pre-smoothing steps to be performed
     169            0 :                 void set_num_presmooth(int num) {m_numPreSmooth = num;}
     170              : 
     171              :         ///     sets the number of post-smoothing steps to be performed
     172            0 :                 void set_num_postsmooth(int num) {m_numPostSmooth = num;}
     173              : 
     174              :         ///     sets the smoother that is used
     175            0 :                 void set_smoother(SmartPtr<ILinearIterator<vector_type> > smoother)
     176            0 :                         {set_presmoother(smoother); set_postsmoother(smoother);}
     177              : 
     178              :         ///     sets the pre-smoother that is used
     179            0 :                 void set_presmoother(SmartPtr<ILinearIterator<vector_type> > smoother)
     180            0 :                         {m_spPreSmootherPrototype = smoother;}
     181              : 
     182              :         ///     sets the post-smoother that is used
     183            0 :                 void set_postsmoother(SmartPtr<ILinearIterator<vector_type> > smoother)
     184            0 :                         {m_spPostSmootherPrototype = smoother;}
     185              : 
     186              :         ///     sets the transfer operator
     187            0 :                 void set_transfer(SmartPtr<ITransferOperator<TDomain, TAlgebra> > P)
     188            0 :                         {set_prolongation(P); set_restriction(P);}
     189              : 
     190              :         ///     sets the prolongation operator
     191            0 :                 void set_prolongation(SmartPtr<ITransferOperator<TDomain, TAlgebra> > P)
     192            0 :                         {m_spProlongationPrototype = P;}
     193              : 
     194              :         ///     sets the restriction operator
     195            0 :                 void set_restriction(SmartPtr<ITransferOperator<TDomain, TAlgebra> > P)
     196            0 :                         {m_spRestrictionPrototype = P;}
     197              : 
     198              :         ///     sets the projection operator
     199            0 :                 void set_projection(SmartPtr<ITransferOperator<TDomain, TAlgebra> > P)
     200            0 :                         {m_spProjectionPrototype = P;}
     201              : 
     202              :         ///     clears all transfer post process
     203            0 :                 void clear_transfer_post_process()
     204            0 :                         {m_vspProlongationPostProcess.clear(); m_vspRestrictionPostProcess.clear();}
     205              : 
     206              :         ///     add prolongation post process
     207            0 :                 void add_prolongation_post_process(SmartPtr<ITransferPostProcess<TDomain, TAlgebra> > PP)
     208            0 :                         {m_vspProlongationPostProcess.push_back(PP);}
     209              : 
     210              :         ///     add restriction post process
     211            0 :                 void add_restriction_post_process(SmartPtr<ITransferPostProcess<TDomain, TAlgebra> > PP)
     212            0 :                         {m_vspRestrictionPostProcess.push_back(PP);}
     213              : 
     214              :         ///////////////////////////////////////////////////////////////////////////
     215              :         //      Linear Solver interface methods
     216              :         ///////////////////////////////////////////////////////////////////////////
     217              : 
     218              :         ///     name
     219            0 :                 virtual const char* name() const {return "Geometric MultiGrid";}
     220              : 
     221              :         ///     returns if parallel solving is supported
     222            0 :                 virtual bool supports_parallel() const
     223              :                 {
     224            0 :                         if(!m_spPreSmootherPrototype->supports_parallel())
     225              :                                 return false;
     226            0 :                         if(!m_spPostSmootherPrototype->supports_parallel())
     227              :                                 return false;
     228              :                         return true;
     229              :                 }
     230              : 
     231              :         ///     returns information about configuration parameters
     232              :                 virtual std::string config_string() const;
     233              : 
     234              :         /// Prepare for Operator J(u) and linearization point u (current solution)
     235              :                 virtual bool init(SmartPtr<ILinearOperator<vector_type> > J, const vector_type& u);
     236              : 
     237              :         ///     Prepare for Linear Operator L
     238              :                 virtual bool init(SmartPtr<ILinearOperator<vector_type> > L);
     239              : 
     240              :         ///     does not call init on base-solver during initialization
     241              :         /**     Use this method with care. It can be useful e.g. during repeated
     242              :          * adaption of a linear problem, where the base-level does not change.
     243              :          * \warning     Currently only works in serial environments
     244              :          * \note        This only affects calls to init(...) and has no effect on
     245              :          *                      direct calls to init_base_solver.
     246              :          * \{ */
     247              :                 virtual void ignore_init_for_base_solver(bool ignore);
     248              :                 virtual bool ignore_init_for_base_solver() const;
     249              :         /** \\ */
     250              : 
     251              :         /// If the nonzero entry pattern of the matrix does not change
     252              :         /// (which is usually the case in time-dependent problems without space adaptation),
     253              :         /// then keeping the entry structure between re-inits saves time during level matrix assembling.
     254            0 :                 void set_matrix_structure_is_const(bool b) {m_bMatrixStructureIsConst = b;}
     255              : 
     256              :         /// reinit transfer operators
     257              :                 void force_reinit();
     258              : 
     259              :         ///     Compute new correction c = B*d
     260              :                 virtual bool apply(vector_type& c, const vector_type& d);
     261              : 
     262              :         ///     Compute new correction c = B*d and return new defect d := d - A*c
     263              :                 virtual bool apply_update_defect(vector_type& c, vector_type& d);
     264              : 
     265              :         ///     Clone
     266              :                 SmartPtr<ILinearIterator<vector_type> > clone();
     267              : 
     268              :         ///     Destructor
     269              :                 ~AssembledMultiGridCycle();
     270              : 
     271              :         protected:
     272              :         /// compute correction on level and update defect
     273              :                 void lmgc(int lev, int cycleType);
     274              : 
     275              :         ////////////////////////////////////////////////////////////////
     276              :         //      The methods in this section rely on each other and should be called in sequence
     277              :         ///     performs presmoothing on the given level
     278              :                 void presmooth_and_restriction(int lev);
     279              : 
     280              :         ///     performs prolongation to the level above
     281              :                 void prolongation_and_postsmooth(int lev);
     282              : 
     283              :         ///     compute base solver
     284              :                 void base_solve(int lev);
     285              :         //      end of section
     286              :         ////////////////////////////////////////////////////////////////
     287              : 
     288              :         ///     allocates the memory
     289              :                 void init_level_memory(int baseLev, int topLev);
     290              : 
     291              :         ///     initializes common part
     292              :                 void init();
     293              : 
     294              :         ///     initializes the smoother and base solver
     295              :                 void init_smoother();
     296              : 
     297              :         ///     initializes the coarse grid matrices
     298              :                 void assemble_level_operator();
     299              :                 void init_rap_operator();
     300              : 
     301              :         ///     initializes the smoother and base solver
     302              :                 void init_base_solver();
     303              : 
     304              :         ///     initializes the transfers
     305              :                 void init_transfer();
     306              : 
     307              :         ///     initializes the prolongation
     308              :                 void init_projection();
     309              : 
     310              :         ///     assembles the missing matrix part on the coarse level, that must be
     311              :         ///     added if the correction has been computed to ensure a correctly updated
     312              :         ///     defect. (i.e. assembles A^c, with d^f -= A^c * c^c)
     313              :                 void assemble_rim_cpl(const vector_type* u);
     314              :                 void init_rap_rim_cpl();
     315              : 
     316              :         protected:
     317              :         /// operator to invert (surface grid)
     318              :                 ConstSmartPtr<matrix_type> m_spSurfaceMat;
     319              : 
     320              :         ///     Solution on surface grid
     321              :                 const vector_type* m_pSurfaceSol;
     322              : 
     323              :         ///     assembling routine for coarse grid matrices
     324              :                 SmartPtr<IAssemble<TAlgebra> > m_spAss;
     325              : 
     326              :         ///     approximation space for level and surface grid
     327              :                 SmartPtr<ApproximationSpace<TDomain> > m_spApproxSpace;
     328              : 
     329              :         ///     top level (i.e. highest level in hierarchy. This is the surface level
     330              :         ///     in case of non-adaptive refinement)
     331              :                 int m_topLev;
     332              :                 int m_surfaceLev;
     333              : 
     334              :         ///     base level (where exact inverse is computed)
     335              :                 int m_baseLev;
     336              : 
     337              :         ///     cylce type (1 = V-cycle, 2 = W-cylcle, ...)
     338              :                 int m_cycleType;
     339              :                 static const int _V_ = 1;
     340              :                 static const int _W_ = 2;
     341              :                 static const int _F_ = -1;
     342              : 
     343              :         ///     number of Presmooth steps
     344              :                 int m_numPreSmooth;
     345              : 
     346              :         ///     number of Postsmooth steps
     347              :                 int m_numPostSmooth;
     348              : 
     349              :         ///     lowest level containing surface geom obj (proc-locally)
     350              :                 int m_LocalFullRefLevel;
     351              : 
     352              :         ///     grid-view for level vectors
     353              :                 GridLevel::ViewType m_GridLevelType;
     354              : 
     355              :         ///     using RAP-Product (assemble coarse-grid matrices otherwise)
     356              :                 bool m_bUseRAP;
     357              : 
     358              :         ///     flag if smoothing on surface rim
     359              :                 bool m_bSmoothOnSurfaceRim;
     360              : 
     361              :         ///     flag if overlapping communication and computation
     362              :                 bool m_bCommCompOverlap;
     363              : 
     364              :         ///     approximation space revision of cached values
     365              :                 RevisionCounter m_ApproxSpaceRevision;
     366              : 
     367              :         ///     prototype for pre-smoother
     368              :                 SmartPtr<ILinearIterator<vector_type> > m_spPreSmootherPrototype;
     369              : 
     370              :         ///     prototype for post-smoother
     371              :                 SmartPtr<ILinearIterator<vector_type> > m_spPostSmootherPrototype;
     372              : 
     373              :         ///     prototype for projection operator
     374              :                 SmartPtr<ITransferOperator<TDomain, TAlgebra> > m_spProjectionPrototype;
     375              : 
     376              :         ///     prototype for prolongation operator
     377              :                 SmartPtr<ITransferOperator<TDomain, TAlgebra> > m_spProlongationPrototype;
     378              : 
     379              :         ///     prototype for restriction operator
     380              :                 SmartPtr<ITransferOperator<TDomain, TAlgebra> > m_spRestrictionPrototype;
     381              : 
     382              :         ///     prototpe for transfer post process
     383              :                 std::vector<SmartPtr<ITransferPostProcess<TDomain, TAlgebra> > > m_vspProlongationPostProcess;
     384              :                 std::vector<SmartPtr<ITransferPostProcess<TDomain, TAlgebra> > > m_vspRestrictionPostProcess;
     385              : 
     386              :         ///     base solver for the coarse problem
     387              :                 SmartPtr<ILinearOperatorInverse<vector_type> > m_spBaseSolver;
     388              : 
     389              :                 ////////////////////////////////////
     390              :                 // Storage for each grid level
     391              :                 ////////////////////////////////////
     392              : 
     393              :         ///     Structure used to realize Surface to Level mapping
     394              :         /// \{
     395              :                 struct LevelIndex{
     396            0 :                         LevelIndex() : index(-1), indexLower(-1), level(-1), levelLower(-1) {}
     397              :                         size_t index, indexLower;
     398              :                         int level, levelLower;
     399              :                 };
     400              :                 std::vector<LevelIndex> m_vSurfToLevelMap;
     401              : 
     402              :                 template <typename TElem>
     403              :                 void init_index_mappings();
     404              :                 void init_index_mappings();
     405              :         /// \}
     406              : 
     407              :         ///     Structure used to realize Surface to Level mapping
     408              :         /// \{
     409              :                 struct SurfLevelMap{
     410              :                         SurfLevelMap() : levIndex(-1), surfIndex(-1) {}
     411            0 :                         SurfLevelMap(size_t levIndex_, size_t surfIndex_)
     412            0 :                                 : levIndex(levIndex_), surfIndex(surfIndex_) {}
     413              :                         size_t levIndex, surfIndex;
     414              :                 };
     415              :         /// \}
     416              : 
     417              :                 struct LevData
     418              :                 {
     419              :                 ///     Level matrix operator
     420              :                         SmartPtr<MatrixOperator<matrix_type, vector_type> > A;
     421              : 
     422              :                 ///     Smoother
     423              :                         SmartPtr<ILinearIterator<vector_type> > PreSmoother;
     424              :                         SmartPtr<ILinearIterator<vector_type> > PostSmoother;
     425              : 
     426              :                 ///     Projection operator
     427              :                         SmartPtr<ITransferOperator<TDomain, TAlgebra> > Projection;
     428              : 
     429              :                 ///     Transfer operator
     430              :                         SmartPtr<ITransferOperator<TDomain, TAlgebra> > Prolongation;
     431              :                         SmartPtr<ITransferOperator<TDomain, TAlgebra> > Restriction;
     432              : 
     433              :                 ///     vectors needed (sx = no-ghosts [for smoothing], t = for transfer)
     434              :                         SmartPtr<GF> sc, sd, st, t;
     435              : 
     436              :                 ///     maps global indices (including ghosts) to patch indices (no ghosts included).
     437              :                         std::vector<size_t> vMapPatchToGlobal;
     438              : 
     439              :                 /// list of shadowing indices
     440              :                         std::vector<size_t> vShadowing;
     441              : 
     442              :                 /// list of corresponding surface index to shadowing indices
     443              :                         std::vector<size_t> vSurfShadowing;
     444              : 
     445              :                 ///     map surface to level
     446              :                         std::vector<SurfLevelMap> vSurfLevelMap;
     447              : 
     448              :                 ///     missing coarse grid correction
     449              :                         matrix_type RimCpl_Fine_Coarse;
     450              : 
     451              :                 ///     missing coarse grid correction
     452              :                         matrix_type RimCpl_Coarse_Fine;
     453              :                         
     454              :                 /// debugging output information (number of calls of the pre-, postsmoothers, base solver etc)
     455              :                         int n_pre_calls, n_post_calls, n_base_calls, n_restr_calls, n_prolong_calls;
     456              :                 };
     457              : 
     458              :         ///     storage for all level
     459              :                 std::vector<SmartPtr<LevData> > m_vLevData;
     460              : 
     461              :         ///     flag, if to solve base problem in parallel when gathered and (!) parallel possible
     462              :                 bool m_bGatheredBaseIfAmbiguous;
     463              : 
     464              :         ///     flag if using parallel base solver
     465              :                 bool m_bGatheredBaseUsed;
     466              : 
     467              :         ///     flag indicating whether the base-solver shall be initialized during a call to init()
     468              :                 bool m_ignoreInitForBaseSolver;
     469              : 
     470              :         /// flag indicating whether the matrix structure is supposed to stay the same in the next init phase
     471              :                 bool m_bMatrixStructureIsConst;
     472              : 
     473              :         ///     Matrix for gathered base solver
     474              :                 SmartPtr<MatrixOperator<matrix_type, vector_type> > spGatheredBaseMat;
     475              : 
     476              :         ///     vector for gathered base solver
     477              :                 SmartPtr<GF> spGatheredBaseCorr;
     478              : 
     479              :         ///     returns if gathered base master
     480              :                 bool gathered_base_master() const;
     481              : 
     482              :         ///     current surface correction
     483              :                 GF* m_pC;
     484              : 
     485              :         ///     init mapping from noghost -> w/ ghost
     486              :         /// \{
     487              :                 template <typename TElem>
     488              :                 void init_noghost_to_ghost_mapping(     std::vector<size_t>& vNoGhostToGhostMap,
     489              :                                                                                         ConstSmartPtr<DoFDistribution> spNoGhostDD,
     490              :                                                                                         ConstSmartPtr<DoFDistribution> spGhostDD);
     491              :                 void init_noghost_to_ghost_mapping(int lev);
     492              :         /// \}
     493              : 
     494              :         ///     copies vector to smoothing patch using cached mapping
     495              :                 void copy_ghost_to_noghost(SmartPtr<GF> spVecTo,
     496              :                                            ConstSmartPtr<GF> spVecFrom,
     497              :                                            const std::vector<size_t>& vMapPatchToGlobal);
     498              : 
     499              :         ///     copies vector from smoothing patch using cached mapping
     500              :                 void copy_noghost_to_ghost(SmartPtr<GF> spVecTo,
     501              :                                                                    ConstSmartPtr<GF> spVecFrom,
     502              :                                                                    const std::vector<size_t>& vMapPatchToGlobal);
     503              : 
     504              :         ///     copies matrix from smoothing patch using cached mapping
     505              :                 void copy_noghost_to_ghost(SmartPtr<matrix_type> spMatTo,
     506              :                                            ConstSmartPtr<matrix_type> spMatFrom,
     507              :                                            const std::vector<size_t>& vMapPatchToGlobal);
     508              : 
     509              : #ifdef UG_PARALLEL
     510              :         /// communicator
     511              :                 pcl::InterfaceCommunicator<IndexLayout> m_Com;
     512              : #endif
     513              : 
     514              :         public:
     515              :         ///     set debug output
     516              :         /**
     517              :          * If a DebugWriter is passed by this method, the multi grid cycle writes
     518              :          * the level/surface vectors and matrices for debug purposes.
     519              :          *
     520              :          * \param[in]   debugWriter             Debug Writer to use
     521              :          */
     522            0 :                 void set_debug(SmartPtr<GridFunctionDebugWriter<TDomain, TAlgebra> > spDebugWriter)
     523              :                 {
     524            0 :                         m_spDebugWriter = spDebugWriter;
     525            0 :                 }
     526              : 
     527              :         protected:
     528              :         ///     writes debug output for a level vector only on smooth path
     529              :         /**
     530              :          * This method writes the level vector to a debug file, if a debug writer
     531              :          * has been set.
     532              :          *
     533              :          * \param[in]           spGF            Level Vector to write for debug purpose
     534              :          * \param[in]           name            Filename
     535              :          */
     536              :                 inline void write_debug(ConstSmartPtr<GF> spGF, std::string name, int cycleNo = -1);
     537              :                 void write_debug(const GF& rGF, std::string name, int cycleNo = -1);
     538              : 
     539              :         ///     writes debug output for a level matrix only on smooth path
     540              :         /**
     541              :          * This method writes the level matrix to a debug file, if a debug writer
     542              :          * has been set.
     543              :          *
     544              :          * \param[in]           mat                     Level Matrix to write for debug purpose
     545              :          * \param[in]           name            Filename
     546              :          */
     547              :         /// \{
     548              :                 void write_debug(const matrix_type& mat, std::string name,
     549              :                                  const GridLevel& glTo, const GridLevel& glFrom);
     550              :                 void write_debug(const matrix_type& mat, std::string name,
     551              :                                  const GF& rTo, const GF& rFrom);
     552              :         /// \}
     553              :         
     554              :         /// enters a new debugger section for smoothers, base solver etc
     555              :                 void enter_debug_writer_section(GridLevel& orig_gl, const char * sec_name, int lev, int cycleNo = -1, int callNo = -1);
     556              :                 
     557              :         /// leaves the current debugger section
     558              :                 void leave_debug_writer_section(GridLevel& orig_gl);
     559              : 
     560              :         ///     logs a level-data-struct to the terminal
     561              :                 void log_debug_data(int lvl, int cycleNo, std::string name);
     562              : 
     563              :         ///     Calls MGStats::set_defect (if available) with the given parameters
     564              :                 void mg_stats_defect(GF& gf, int lvl, typename mg_stats_type::Stage stage);
     565              : 
     566              :         ///     Debug Writer
     567              :                 SmartPtr<GridFunctionDebugWriter<TDomain, TAlgebra> > m_spDebugWriter;
     568              : 
     569              :         ///     counter for debug, to distinguish the iterations
     570              :                 int m_dbgIterCnt;
     571              : 
     572              :         ///     MGStats are used for debugging to record statistics on individual iterations
     573              :                 SmartPtr<mg_stats_type>   m_mgstats;
     574              : };
     575              : 
     576              : ////////////////////////////////////////////////////////////////////////////////
     577              : // Selections
     578              : ////////////////////////////////////////////////////////////////////////////////
     579              : 
     580              : /// selects all non-shadows, that are adjacent to a shadow on a grid levels
     581              : void SelectNonShadowsAdjacentToShadowsOnLevel(BoolMarker& sel,
     582              :                                                                                    const SurfaceView& surfView,
     583              :                                                                                    int level);
     584              : 
     585              : } // end namespace ug
     586              : 
     587              : // include implementation
     588              : #include "mg_solver_impl.hpp"
     589              : 
     590              : #endif /* __H__UG__LIB_DISC__MULTI_GRID_SOLVER__MG_SOLVER__ */
        

Generated by: LCOV version 2.0-1