Line data Source code
1 : /*
2 : * Copyright (c) 2012-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__FUNCTION_SPACS__LOCAL_TRANSFER_INTERFACE__
34 : #define __H__UG__LIB_DISC__FUNCTION_SPACS__LOCAL_TRANSFER_INTERFACE__
35 :
36 : #include "lib_disc/domain.h"
37 : #include "lib_disc/local_finite_element/local_finite_element_id.h"
38 :
39 : namespace ug{
40 :
41 : ////////////////////////////////////////////////////////////////////////////////
42 : // Value Access Base class
43 : ////////////////////////////////////////////////////////////////////////////////
44 :
45 : class TransferValueAccessor
46 : {
47 : public:
48 : virtual void access_inner(GridObject* elem) = 0;
49 :
50 : virtual void access_closure(GridObject* elem) = 0;
51 :
52 : virtual ~TransferValueAccessor() {}
53 :
54 : const number& operator[](size_t i) const {
55 : UG_ASSERT(i < m_Val.size(), "Wrong index "<<i<<" (size: "<<m_Val.size()<<")");
56 : return *m_Val[i];
57 : }
58 :
59 : number& operator[](size_t i) {
60 : UG_ASSERT(i < m_Val.size(), "Wrong index "<<i<<" (size: "<<m_Val.size()<<")");
61 0 : return *m_Val[i];
62 : }
63 :
64 : size_t size() const {return m_Val.size();}
65 :
66 : protected:
67 : std::vector<number*> m_Val;
68 : };
69 :
70 : ////////////////////////////////////////////////////////////////////////////////
71 : // Element Prolongation
72 : ////////////////////////////////////////////////////////////////////////////////
73 :
74 : template <typename TDomain>
75 : class IElemProlongation
76 : {
77 : public:
78 0 : virtual void init(ConstSmartPtr<TDomain> spDomain,
79 : SmartPtr<TransferValueAccessor> vValueChild,
80 : SmartPtr<TransferValueAccessor> vValueParent)
81 : {
82 0 : m_spDomain = spDomain;
83 0 : m_spGrid = m_spDomain->grid();
84 0 : m_vValueChild = vValueChild;
85 0 : m_vValueParent = vValueParent;
86 0 : }
87 :
88 : virtual bool perform_prolongation_on(GridBaseObjectId gbo) = 0;
89 :
90 : virtual void prolongate(Vertex* parent) = 0;
91 : virtual void prolongate(Edge* parent) = 0;
92 : virtual void prolongate(Face* parent) = 0;
93 : virtual void prolongate(Volume* parent) = 0;
94 :
95 0 : virtual ~IElemProlongation() {}
96 :
97 : protected:
98 : ConstSmartPtr<TDomain> m_spDomain;
99 : ConstSmartPtr<MultiGrid> m_spGrid;
100 :
101 : SmartPtr<TransferValueAccessor> m_vValueChild;
102 : SmartPtr<TransferValueAccessor> m_vValueParent;
103 : };
104 :
105 :
106 : template <typename TDomain, typename TImpl>
107 : class ElemProlongationBase : public IElemProlongation<TDomain>
108 : {
109 : public:
110 0 : virtual void prolongate(Vertex* parent){
111 0 : this->getImpl().prolongate(parent, *this->m_vValueChild, *this->m_vValueParent);
112 0 : }
113 0 : virtual void prolongate(Edge* parent){
114 0 : this->getImpl().prolongate(parent, *this->m_vValueChild, *this->m_vValueParent);
115 0 : }
116 0 : virtual void prolongate(Face* parent){
117 0 : this->getImpl().prolongate(parent, *this->m_vValueChild, *this->m_vValueParent);
118 0 : }
119 0 : virtual void prolongate(Volume* parent){
120 0 : this->getImpl().prolongate(parent, *this->m_vValueChild, *this->m_vValueParent);
121 0 : }
122 :
123 : protected:
124 : /// access to implementation
125 : TImpl& getImpl() {return static_cast<TImpl&>(*this);}
126 :
127 : /// const access to implementation
128 : const TImpl& getImpl() const {return static_cast<const TImpl&>(*this);}
129 : };
130 :
131 :
132 : template <typename TDomain>
133 : SmartPtr<IElemProlongation<TDomain> >
134 : GetStandardElementProlongation(const LFEID& lfeid);
135 :
136 :
137 : ////////////////////////////////////////////////////////////////////////////////
138 : // Element Restriction
139 : ////////////////////////////////////////////////////////////////////////////////
140 :
141 : template <typename TDomain>
142 : class IElemRestriction
143 : {
144 : public:
145 0 : virtual void init(ConstSmartPtr<TDomain> spDomain,
146 : SmartPtr<TransferValueAccessor> vValueChild,
147 : SmartPtr<TransferValueAccessor> vValueParent)
148 : {
149 0 : m_spDomain = spDomain;
150 0 : m_spGrid = m_spDomain->grid();
151 0 : m_vValueChild = vValueChild;
152 0 : m_vValueParent = vValueParent;
153 0 : }
154 :
155 : virtual bool perform_restriction_on(GridBaseObjectId gbo) = 0;
156 :
157 : virtual void do_restrict(Vertex* parent) = 0;
158 : virtual void do_restrict(Edge* parent) = 0;
159 : virtual void do_restrict(Face* parent) = 0;
160 : virtual void do_restrict(Volume* parent) = 0;
161 :
162 0 : virtual ~IElemRestriction() {}
163 :
164 : protected:
165 : ConstSmartPtr<TDomain> m_spDomain;
166 : ConstSmartPtr<MultiGrid> m_spGrid;
167 :
168 : SmartPtr<TransferValueAccessor> m_vValueChild;
169 : SmartPtr<TransferValueAccessor> m_vValueParent;
170 : };
171 :
172 :
173 : template <typename TDomain, typename TImpl>
174 : class ElemRestrictionBase : public IElemRestriction<TDomain>
175 : {
176 : public:
177 0 : virtual void do_restrict(Vertex* parent){
178 0 : this->getImpl().do_restrict(parent, *this->m_vValueChild, *this->m_vValueParent);
179 0 : }
180 0 : virtual void do_restrict(Edge* parent){
181 0 : this->getImpl().do_restrict(parent, *this->m_vValueChild, *this->m_vValueParent);
182 0 : }
183 0 : virtual void do_restrict(Face* parent){
184 0 : this->getImpl().do_restrict(parent, *this->m_vValueChild, *this->m_vValueParent);
185 0 : }
186 0 : virtual void do_restrict(Volume* parent){
187 0 : this->getImpl().do_restrict(parent, *this->m_vValueChild, *this->m_vValueParent);
188 0 : }
189 :
190 : protected:
191 : /// access to implementation
192 0 : TImpl& getImpl() {return static_cast<TImpl&>(*this);}
193 :
194 : /// const access to implementation
195 : const TImpl& getImpl() const {return static_cast<const TImpl&>(*this);}
196 : };
197 :
198 :
199 : template <typename TDomain>
200 : SmartPtr<IElemRestriction<TDomain> >
201 : GetStandardElementRestriction(const LFEID& lfeid);
202 :
203 :
204 : } // end namespace ug
205 :
206 : #endif /* __H__UG__LIB_DISC__FUNCTION_SPACS__LOCAL_TRANSFER_INTERFACE__ */
|