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__OPERATOR__LINEAR_OPERATOR__STD_TRANSFER__
34 : #define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__STD_TRANSFER__
35 :
36 : // extern headers
37 : #include <iostream>
38 :
39 : // other ug4 modules
40 : #include "common/common.h"
41 : #include "transfer_interface.h"
42 : #include "lib_algebra/operator/debug_writer.h"
43 :
44 : #ifdef UG_PARALLEL
45 : #include "lib_disc/parallelization/parallelization_util.h"
46 : #endif
47 :
48 : namespace ug{
49 :
50 : /// Standard Prolongation Operator
51 : /** By default a special optimization is performed for p1-lagrange-elements.
52 : * This optimization is only valid if all elements have been refined with
53 : * standard refinement rules. If closure elements are generated, this optimization
54 : * has to be deactivated (use StdTransfer::enable_p1_lagrange_optimization(false)).
55 : */
56 : template <typename TDomain, typename TAlgebra>
57 : class StdTransfer :
58 : virtual public ITransferOperator<TDomain, TAlgebra>
59 : {
60 : public:
61 : /// Type of base class
62 : typedef ITransferOperator<TDomain, TAlgebra> base_type;
63 :
64 : /// Type of Algebra
65 : typedef TAlgebra algebra_type;
66 :
67 : /// Type of Vector
68 : typedef typename TAlgebra::vector_type vector_type;
69 :
70 : /// Type of Matrix
71 : typedef typename TAlgebra::matrix_type matrix_type;
72 :
73 : /// Type of Domain
74 : typedef TDomain domain_type;
75 :
76 : /// Type of GridFunction
77 : typedef GridFunction<TDomain, TAlgebra> GF;
78 :
79 : public:
80 : /// Default constructor
81 0 : StdTransfer() : ITransferOperator<TDomain, TAlgebra>(),
82 0 : m_p1LagrangeOptimizationEnabled(true),
83 0 : m_dampRes(1.0), m_dampProl(1.0),
84 0 : bCached(true), m_bUseTransposed(true),
85 0 : m_spDebugWriter(NULL)
86 : {};
87 :
88 : /// virtual destructor
89 0 : virtual ~StdTransfer(){};
90 :
91 : /// set restriction damping (only applied on vector operation, not (!!) in assembled matrices)
92 0 : void set_restriction_damping(number damp) {m_dampRes = damp;}
93 :
94 : /// set prolongation damping (only applied on vector operation, not (!!) in assembled matrices)
95 0 : void set_prolongation_damping(number damp) {m_dampProl = damp;}
96 :
97 : /// set debug writer
98 0 : void set_debug(SmartPtr<IDebugWriter<TAlgebra> > spDebugWriter) {
99 0 : m_spDebugWriter = spDebugWriter;
100 0 : }
101 :
102 : /// enables/disables an assembling optimization for p1-lagrange elements
103 : /** The optimization is enabled by default. It can however only be used,
104 : * if all elements are refined with their standard refinement rule. If one
105 : * uses anisotropic refinement or refinement with closure, the optimization
106 : * should be disabled.
107 : * \todo The normal assembling strategy should be optimized in such a way
108 : * that the p1-lagrange optimization is no longer required. This
109 : * however involves something like a ref-type-hash for each element,
110 : * which returns a unique number based on the types and order of children.*/
111 0 : void enable_p1_lagrange_optimization(bool enable) {m_p1LagrangeOptimizationEnabled = enable;}
112 0 : bool p1_lagrange_optimization_enabled() const {return m_p1LagrangeOptimizationEnabled;}
113 :
114 : /// sets if restriction and prolongation are transposed
115 0 : void set_use_transposed(bool bTransposed) {m_bUseTransposed = bTransposed;}
116 :
117 : public:
118 : /// Set levels
119 0 : virtual void set_levels(GridLevel coarseLevel, GridLevel fineLevel) {}
120 :
121 : /// initialize the operator
122 0 : virtual void init() {}
123 :
124 : /// returns new instance with same setting
125 : virtual SmartPtr<ITransferOperator<TDomain, TAlgebra> > clone();
126 :
127 : /// apply Operator, interpolate function
128 0 : virtual void prolongate(vector_type& uFine, const vector_type& uCoarse){
129 0 : GF* pFine = dynamic_cast<GF*>(&uFine);
130 0 : const GF* pCoarse = dynamic_cast<const GF*>(&uCoarse);
131 0 : if(!pFine || !pCoarse)
132 0 : UG_THROW("StdTransfer: fine and coarse vectors expected to be "
133 : "a grid function.");
134 0 : prolongate(*pFine, *pCoarse);
135 0 : }
136 :
137 : /// apply transposed Operator, restrict function
138 0 : virtual void do_restrict(vector_type& uCoarse, const vector_type& uFine){
139 0 : const GF* pFine = dynamic_cast<const GF*>(&uFine);
140 0 : GF* pCoarse = dynamic_cast<GF*>(&uCoarse);
141 0 : if(!pFine || !pCoarse)
142 0 : UG_THROW("StdTransfer: fine and coarse vectors expected to be "
143 : "a grid function.");
144 0 : do_restrict(*pCoarse, *pFine);
145 0 : }
146 :
147 : public:
148 : /// returns prolongation as a matrix
149 : virtual SmartPtr<matrix_type>
150 : prolongation(const GridLevel& fineGL, const GridLevel& coarseGL,
151 : ConstSmartPtr<ApproximationSpace<TDomain> > spApproxSpace);
152 :
153 : /// returns restriction as a matrix
154 : virtual SmartPtr<matrix_type>
155 : restriction(const GridLevel& coarseGL, const GridLevel& fineGL,
156 : ConstSmartPtr<ApproximationSpace<TDomain> > spApproxSpace);
157 :
158 : /// apply operator to a grid function
159 : void prolongate(GF& uFine, const GF& uCoarse);
160 :
161 : /// apply operator to a grid function
162 : void do_restrict(GF& uCoarse, const GF& uFine);
163 :
164 : protected:
165 : /// debug writing of matrix
166 : void write_debug(const matrix_type& mat, std::string name,
167 : const GridLevel& glTo, const GridLevel& glFrom);
168 :
169 : template <typename TChild>
170 : void assemble_restriction(matrix_type& mat,
171 : const DoFDistribution& coarseDD,
172 : const DoFDistribution& fineDD,
173 : ConstSmartPtr<TDomain> spDomain);
174 : void assemble_restriction(matrix_type& mat,
175 : const DoFDistribution& coarseDD,
176 : const DoFDistribution& fineDD,
177 : ConstSmartPtr<TDomain> spDomain);
178 :
179 : template <typename TChild>
180 : void assemble_prolongation(matrix_type& mat,
181 : const DoFDistribution& fineDD,
182 : const DoFDistribution& coarseDD,
183 : ConstSmartPtr<TDomain> spDomain);
184 : void assemble_prolongation(matrix_type& mat,
185 : const DoFDistribution& fineDD,
186 : const DoFDistribution& coarseDD,
187 : ConstSmartPtr<TDomain> spDomain);
188 :
189 : void assemble_prolongation_p1(matrix_type& mat,
190 : const DoFDistribution& fineDD,
191 : const DoFDistribution& coarseDD);
192 :
193 : protected:
194 : /// struct to distinguish already assembled operators
195 : struct TransferKey{
196 0 : TransferKey(const GridLevel& toGL_, const GridLevel& fromGL_,
197 : const RevisionCounter& revCnt_)
198 0 : : toGL(toGL_), fromGL(fromGL_), revCnt(revCnt_) {}
199 : GridLevel toGL, fromGL;
200 : RevisionCounter revCnt;
201 :
202 0 : bool operator<(const TransferKey& other) const {
203 0 : if(revCnt != other.revCnt) return revCnt < other.revCnt;
204 : if(toGL != other.toGL) return toGL < other.toGL;
205 : return fromGL < other.fromGL;
206 : }
207 : };
208 :
209 : typedef std::map<TransferKey, SmartPtr<matrix_type> > TransferMap;
210 : TransferMap m_mRestriction;
211 : TransferMap m_mProlongation;
212 :
213 0 : void remove_outdated(TransferMap& map, const RevisionCounter& revCnt) {
214 : typedef typename TransferMap::iterator iterator;
215 0 : for(iterator iter = map.begin(); iter != map.end();)
216 : {
217 : const RevisionCounter& cnt = iter->first.revCnt;
218 0 : if((cnt.obj() == revCnt.obj()) && (cnt != revCnt)){
219 : map.erase(iter++);
220 : } else {
221 : ++iter;
222 : }
223 : }
224 0 : }
225 :
226 : protected:
227 : /// list of post processes
228 : using base_type::m_vConstraint;
229 :
230 : /// flag for p1-lagrange-optimization
231 : bool m_p1LagrangeOptimizationEnabled;
232 :
233 : /// damping parameter
234 : number m_dampRes;
235 : number m_dampProl;
236 :
237 : /// flag if cached (matrix) transfer used
238 : bool bCached;
239 :
240 : /// flag if transposed is used
241 : bool m_bUseTransposed;
242 :
243 : /// debug writer
244 : SmartPtr<IDebugWriter<TAlgebra> > m_spDebugWriter;
245 : };
246 :
247 : } // end namespace ug
248 :
249 : #include "std_transfer_impl.h"
250 :
251 : #endif /* __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__STD_TRANSFER__ */
|