Line data Source code
1 : /*
2 : * Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
3 : * Authors: Andreas Vogel, Christian Wehner
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__
34 : #define __H__UG__LIB_DISC__FUNCTION_SPACS__LOCAL_TRANSFER__
35 :
36 : #include "local_transfer_interface.h"
37 : #include "lib_disc/local_finite_element/local_finite_element_provider.h"
38 : #include "lib_disc/domain_util.h"
39 : #include "lib_disc/reference_element/reference_mapping_provider.h"
40 : #include "lib_disc/function_spaces/dof_position_util.h"
41 :
42 : namespace ug{
43 :
44 :
45 : template <typename TDomain>
46 : class PiecewiseConstantElemTransfer
47 : : public ElemProlongationBase<TDomain, PiecewiseConstantElemTransfer<TDomain> >
48 : , public ElemRestrictionBase<TDomain, PiecewiseConstantElemTransfer<TDomain> >
49 : {
50 : public:
51 0 : PiecewiseConstantElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
52 :
53 0 : virtual bool perform_prolongation_on(GridBaseObjectId gbo) {
54 0 : if(m_lfeid.dim() == gbo) return true;
55 : return false;
56 : }
57 :
58 : template <typename TElem>
59 0 : void prolongate(TElem* parent,
60 : TransferValueAccessor& vValueChild,
61 : TransferValueAccessor& vValueParent)
62 : {
63 : const MultiGrid* mg = IElemProlongation<TDomain>::m_spGrid.get();
64 0 : const int numChild = mg->num_children<TElem>(parent);
65 :
66 0 : vValueParent.access_inner(parent);
67 0 : for(int c = 0; c < numChild; ++c){
68 0 : TElem* child = mg->get_child<TElem>(parent, c);
69 :
70 0 : vValueChild.access_inner(child);
71 0 : vValueChild[0] = vValueParent[0];
72 : }
73 0 : }
74 :
75 0 : virtual bool perform_restriction_on(GridBaseObjectId gbo) {
76 0 : if(m_lfeid.dim() == gbo) return true;
77 : return false;
78 : }
79 :
80 : template <typename TElem>
81 0 : void do_restrict(TElem* parent,
82 : TransferValueAccessor& vValueChild,
83 : TransferValueAccessor& vValueParent)
84 : {
85 : const MultiGrid* mg = IElemRestriction<TDomain>::m_spGrid.get();
86 0 : const int numChild = mg->num_children<TElem>(parent);
87 :
88 0 : vValueParent.access_inner(parent);
89 : UG_ASSERT(vValueParent.size()==1, "Exactly one DoF for Piecewise-Constant")
90 0 : vValueParent[0] = 0.0;
91 :
92 0 : for(int c = 0; c < numChild; ++c){
93 0 : TElem* child = mg->get_child<TElem>(parent, c);
94 :
95 0 : vValueChild.access_inner(child);
96 : UG_ASSERT(vValueChild.size()==1, "Exactly one DoF for Piecewise-Constant")
97 0 : vValueParent[0] += vValueChild[0];
98 : }
99 0 : vValueParent[0] *= (1./numChild);
100 0 : }
101 :
102 : protected:
103 : LFEID m_lfeid;
104 : };
105 :
106 : template <typename TDomain>
107 : class P1LagrangeElemTransfer
108 : : public ElemProlongationBase<TDomain, P1LagrangeElemTransfer<TDomain> >
109 : , public ElemRestrictionBase<TDomain, P1LagrangeElemTransfer<TDomain> >
110 : {
111 : public:
112 0 : P1LagrangeElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
113 :
114 0 : virtual bool perform_prolongation_on(GridBaseObjectId gbo) {
115 0 : if(m_lfeid.dim() < gbo) return false;
116 : return true;
117 : }
118 :
119 : template <typename TElem>
120 0 : void prolongate(TElem* parent,
121 : TransferValueAccessor& vValueChild,
122 : TransferValueAccessor& vValueParent)
123 : {
124 : const MultiGrid* mg = IElemProlongation<TDomain>::m_spGrid.get();
125 : const int numChild = mg->num_children<Vertex>(parent);
126 : if(numChild == 0) return;
127 : if(numChild != 1) UG_THROW("Max child Vertex must be 1");
128 :
129 : Vertex* child = mg->get_child<Vertex>(parent, 0);
130 :
131 0 : vValueChild.access_inner(child);
132 0 : vValueParent.access_closure(parent);
133 :
134 0 : if (vValueChild.size() == 0)
135 : return; // current fct not defined on child
136 :
137 0 : const int numVrt = vValueParent.size();
138 : UG_ASSERT(numVrt > 0, "No function value found on parent vertices.");
139 :
140 0 : vValueChild[0] = vValueParent[0];
141 0 : for(int i = 1; i < numVrt; ++i)
142 0 : vValueChild[0] += vValueParent[i];
143 0 : vValueChild[0] *= 1.0/numVrt;
144 : }
145 :
146 0 : virtual bool perform_restriction_on(GridBaseObjectId gbo){
147 0 : if(gbo == VERTEX) return true;
148 : return false;
149 : }
150 :
151 : // the following line silences -Woverloaded-virtual
152 : using ElemRestrictionBase<TDomain, P1LagrangeElemTransfer<TDomain> >::do_restrict;
153 0 : void do_restrict(GridObject* parent,
154 : TransferValueAccessor& vValueChild,
155 : TransferValueAccessor& vValueParent)
156 : {
157 0 : UG_THROW("Should not be called.")
158 : }
159 :
160 0 : void do_restrict(Vertex* parent,
161 : TransferValueAccessor& vValueChild,
162 : TransferValueAccessor& vValueParent)
163 : {
164 : const MultiGrid* mg = IElemRestriction<TDomain>::m_spGrid.get();
165 : const int numChild = mg->num_children<Vertex>(parent);
166 0 : if(numChild != 1) UG_THROW("Num child Vertex must be 1");
167 :
168 : Vertex* child = mg->get_child<Vertex>(parent, 0);
169 :
170 0 : vValueChild.access_inner(child);
171 0 : vValueParent.access_closure(parent);
172 :
173 0 : if (vValueParent.size() == 0)
174 : return; // current fct not defined on parent
175 :
176 : UG_ASSERT(vValueChild.size() > 0, "No function value found on child vertex.");
177 :
178 0 : vValueParent[0] = vValueChild[0];
179 : }
180 :
181 : protected:
182 : LFEID m_lfeid;
183 : };
184 :
185 :
186 : template <typename TDomain>
187 : class StdLagrangeElemTransfer
188 : : public ElemProlongationBase<TDomain,StdLagrangeElemTransfer<TDomain> >
189 : , public ElemRestrictionBase<TDomain,StdLagrangeElemTransfer<TDomain> >
190 : {
191 : public:
192 : /// world dimension
193 : static const int dim = TDomain::dim;
194 :
195 : public:
196 0 : StdLagrangeElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
197 :
198 0 : virtual bool perform_prolongation_on(GridBaseObjectId gbo) {
199 0 : if(m_lfeid.order() == 1 && gbo != VERTEX) return false;
200 0 : if(m_lfeid.dim() < gbo) return false;
201 : return true;
202 : }
203 :
204 :
205 : template <typename TParent, typename TChild>
206 0 : void prolongate(TParent* parent,
207 : TransferValueAccessor& vValueChild,
208 : TransferValueAccessor& vValueParent,
209 : const LocalShapeFunctionSet<TParent::dim>& lsfs,
210 : const std::vector<MathVector<dim> >& vCornerParent,
211 : const ReferenceObjectID parentRoid)
212 : {
213 : const MultiGrid* mg = IElemProlongation<TDomain>::m_spGrid.get();
214 : static const int pdim = TParent::dim;
215 :
216 : // number of children of the base type
217 0 : const int numChild = mg->num_children<TChild>(parent);
218 :
219 : // loop all children
220 0 : for(int c = 0; c < numChild; ++c)
221 : {
222 : // get child
223 0 : TChild* child = mg->get_child<TChild>(parent, c);
224 :
225 : // access the values
226 0 : vValueChild.access_inner(child);
227 :
228 : // global positions of fine dofs
229 : std::vector<MathVector<dim> > vDoFPos;
230 0 : InnerDoFPosition(vDoFPos, child, *IElemProlongation<TDomain>::m_spDomain, m_lfeid);
231 :
232 : // get Reference Mapping
233 : DimReferenceMapping<pdim, dim>& map =
234 : ReferenceMappingProvider::get<pdim, dim>(parentRoid, vCornerParent);
235 :
236 : // get local position of DoF
237 0 : std::vector<MathVector<pdim> > vLocPos(vDoFPos.size(), 0.0);
238 0 : map.global_to_local(vLocPos, vDoFPos);
239 :
240 : // evaluate coarse shape fct at fine local point
241 : std::vector<std::vector<number> > vvShape;
242 0 : lsfs.shapes(vvShape, vLocPos);
243 :
244 : // sum up interpolated values
245 0 : for(size_t csh = 0; csh < vValueChild.size(); ++csh){
246 0 : vValueChild[csh] = 0.0;
247 0 : for(size_t psh = 0; psh < vValueParent.size(); ++psh)
248 0 : vValueChild[csh] += vvShape[csh][psh] * vValueParent[psh];
249 : }
250 : }
251 0 : }
252 :
253 :
254 : template <typename TElem>
255 0 : void prolongate(TElem* parent,
256 : TransferValueAccessor& vValueChild,
257 : TransferValueAccessor& vValueParent)
258 : {
259 : static const int pdim = TElem::dim;
260 :
261 : // get reference object ids
262 0 : const ReferenceObjectID parentRoid = parent->reference_object_id();
263 :
264 : // access the values
265 0 : vValueParent.access_closure(parent);
266 :
267 : // get vertices if required by some prolongation
268 : std::vector<MathVector<dim> > vCornerParent;
269 : CollectCornerCoordinates(vCornerParent, *parent, *IElemProlongation<TDomain>::m_spDomain);
270 :
271 : // get local finite element trial spaces
272 : const LocalShapeFunctionSet<pdim>& lsfs =
273 0 : LocalFiniteElementProvider::get<pdim>(parentRoid, m_lfeid);
274 :
275 : // prolongate to all children with dimension pdim and lower
276 : switch(pdim){
277 0 : case 3: prolongate<TElem,Volume>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
278 0 : case 2: prolongate<TElem,Face>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
279 0 : case 1: prolongate<TElem,Edge>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
280 0 : prolongate<TElem,Vertex>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
281 : break;
282 0 : default: UG_THROW("Dimension "<<pdim<<" not supported");
283 : }
284 0 : }
285 :
286 0 : virtual void prolongate(Vertex* parent)
287 : {
288 : const MultiGrid* mg = IElemProlongation<TDomain>::m_spGrid.get();
289 : TransferValueAccessor& vValueChild = *IElemProlongation<TDomain>::m_vValueChild;
290 : TransferValueAccessor& vValueParent = *IElemProlongation<TDomain>::m_vValueParent;
291 :
292 : // access the values
293 0 : vValueParent.access_closure(parent);
294 :
295 : // number of children of the base type
296 : const int numChild = mg->num_children<Vertex>(parent);
297 : if(numChild == 0) return;
298 : if(numChild != 1) UG_THROW("Num child Vertex must be 1");
299 :
300 : // get child
301 : Vertex* child = mg->get_child<Vertex>(parent, 0);
302 :
303 : // access the values
304 0 : vValueChild.access_inner(child);
305 :
306 : // copy value
307 0 : vValueChild[0] = vValueParent[0];
308 : }
309 :
310 :
311 0 : virtual bool perform_restriction_on(GridBaseObjectId gbo) {
312 0 : if(m_lfeid.dim() < gbo) return false;
313 : return true;
314 : }
315 :
316 : template <typename TElem>
317 0 : void do_restrict(TElem* parent,
318 : TransferValueAccessor& vValueChild,
319 : TransferValueAccessor& vValueParent)
320 : {
321 : // access the values
322 0 : vValueParent.access_inner(parent);
323 :
324 : // global positions of fine dofs
325 : std::vector<MathVector<dim> > vDoFPos;
326 0 : InnerDoFPosition(vDoFPos, parent, *IElemRestriction<TDomain>::m_spDomain, m_lfeid);
327 :
328 : const MultiGrid* mg = IElemRestriction<TDomain>::m_spGrid.get();
329 :
330 : // number of children of the base type
331 0 : const int numChild = mg->num_children<TElem>(parent);
332 :
333 : // loop all children
334 0 : for(size_t ip = 0; ip < vDoFPos.size(); ++ip)
335 : {
336 : int c = 0;
337 0 : for(; c < numChild; ++c)
338 : {
339 : // get child
340 0 : TElem* child = mg->get_child<TElem>(parent, c);
341 :
342 : // global positions of fine dofs
343 : std::vector<MathVector<dim> > vChildDoFPos;
344 0 : DoFPosition(vChildDoFPos, child, *IElemRestriction<TDomain>::m_spDomain, m_lfeid);
345 :
346 : // check if pos matches
347 : int cip = -1;
348 0 : for(size_t i = 0; i < vChildDoFPos.size(); ++i)
349 0 : if(VecDistanceSq(vChildDoFPos[i], vDoFPos[ip]) < 1e-8)
350 0 : cip = i;
351 :
352 : // if not found in this child, continue
353 0 : if(cip == -1) continue;
354 :
355 : // access the values
356 0 : vValueChild.access_closure(child);
357 :
358 : // sum up interpolated values
359 0 : vValueParent[ip] = vValueChild[cip];
360 :
361 : break;
362 : }
363 : // check that DoF found
364 0 : if(c == numChild)
365 0 : UG_THROW("Lagrange-Element-Coarsening: Cannot find DoF position "
366 : "in child elems for "<<ip<<"'s Pos: "<<vDoFPos[ip])
367 : }
368 :
369 0 : }
370 :
371 :
372 0 : virtual void do_restrict(Vertex* parent)
373 : {
374 : const MultiGrid* mg = IElemRestriction<TDomain>::m_spGrid.get();
375 : const int numChild = mg->num_children<Vertex>(parent);
376 0 : if(numChild != 1) UG_THROW("Num child Vertex must be 1");
377 :
378 : TransferValueAccessor& vValueChild = *IElemRestriction<TDomain>::m_vValueChild;
379 : TransferValueAccessor& vValueParent = *IElemRestriction<TDomain>::m_vValueParent;
380 :
381 : Vertex* child = mg->get_child<Vertex>(parent, 0);
382 :
383 0 : vValueChild.access_inner(child);
384 0 : vValueParent.access_inner(parent);
385 :
386 0 : vValueParent[0] = vValueChild[0];
387 0 : }
388 :
389 :
390 : protected:
391 : LFEID m_lfeid;
392 : };
393 :
394 :
395 : template <typename TDomain>
396 : class CrouzeixRaviartElemTransfer
397 : : public ElemProlongationBase<TDomain, CrouzeixRaviartElemTransfer<TDomain> >
398 : , public ElemRestrictionBase<TDomain, CrouzeixRaviartElemTransfer<TDomain> >
399 : {
400 : public:
401 : /// world dimension
402 : static const int dim = TDomain::dim;
403 :
404 : public:
405 0 : CrouzeixRaviartElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
406 :
407 0 : virtual bool perform_prolongation_on(GridBaseObjectId gbo) {
408 : // prolongation from elems that have a side as child
409 0 : if(m_lfeid.dim()-1 == gbo) return true;
410 0 : if(m_lfeid.dim() == gbo) return true;
411 : return false;
412 : }
413 :
414 : // the following line silences -Woverloaded-virtual
415 : using ElemProlongationBase<TDomain, CrouzeixRaviartElemTransfer<TDomain> >::prolongate;
416 0 : void prolongate(Vertex* parent,
417 : TransferValueAccessor& vValueChild,
418 : TransferValueAccessor& vValueParent)
419 : {
420 0 : UG_THROW("No implementation for Vertex and Crouzeix-Raviart.")
421 : }
422 :
423 : template <typename TSide>
424 0 : void prolongate(const std::vector<typename TSide::sideof*>& vParentElem,
425 : const std::vector<TSide*>& vChildSide,
426 : TransferValueAccessor& vValueChild,
427 : TransferValueAccessor& vValueParent)
428 : {
429 : // get associated macro elements, this elem is a side of
430 : typedef typename TSide::sideof TElem;
431 0 : if(vParentElem.size() > 2)
432 0 : UG_THROW("More than 2 Elements share a side.")
433 :
434 : // reset values to zero
435 0 : for(size_t c = 0; c < vChildSide.size(); ++c){
436 0 : vValueChild.access_inner(vChildSide[c]);
437 0 : for(size_t csh = 0; csh < vValueChild.size(); ++csh)
438 0 : vValueChild[csh] = 0.0;
439 : }
440 :
441 : // loop parent Elems
442 0 : for(size_t p = 0; p < vParentElem.size(); ++p)
443 : {
444 : // get parent elem
445 0 : TElem* parentElem = vParentElem[p];
446 :
447 : // get reference object ids
448 : static const int pdim = TElem::dim;
449 0 : const ReferenceObjectID parentRoid = parentElem->reference_object_id();
450 :
451 : // access the values
452 0 : vValueParent.access_closure(parentElem);
453 :
454 : // get vertices if required by some prolongation
455 : std::vector<MathVector<dim> > vCornerParent;
456 : CollectCornerCoordinates(vCornerParent, *parentElem, *IElemProlongation<TDomain>::m_spDomain);
457 :
458 : // get local finite element trial spaces
459 : const LocalShapeFunctionSet<pdim>& lsfs =
460 0 : LocalFiniteElementProvider::get<pdim>(parentRoid, m_lfeid);
461 :
462 : // loop child sides
463 0 : for(size_t c = 0; c < vChildSide.size(); ++c)
464 : {
465 : // get child side
466 0 : TSide* childSide = vChildSide[c];
467 :
468 : // access the value
469 0 : vValueChild.access_inner(childSide);
470 :
471 : // global positions of fine dofs
472 : std::vector<MathVector<dim> > vDoFPos;
473 0 : InnerDoFPosition(vDoFPos, childSide, *IElemProlongation<TDomain>::m_spDomain, m_lfeid);
474 :
475 : // get Reference Mapping
476 : DimReferenceMapping<pdim, dim>& map =
477 : ReferenceMappingProvider::get<pdim, dim>(parentRoid, vCornerParent);
478 :
479 : // get local position of DoF
480 0 : std::vector<MathVector<pdim> > vLocPos(vDoFPos.size(), 0.0);
481 0 : map.global_to_local(vLocPos, vDoFPos);
482 :
483 : // evaluate coarse shape fct at fine local point
484 : std::vector<std::vector<number> > vvShape;
485 0 : lsfs.shapes(vvShape, vLocPos);
486 :
487 : // sum up interpolated values
488 0 : for(size_t csh = 0; csh < vValueChild.size(); ++csh){
489 0 : for(size_t psh = 0; psh < vValueParent.size(); ++psh)
490 0 : vValueChild[csh] += (1./vParentElem.size()) *
491 0 : vvShape[csh][psh] * vValueParent[psh];
492 : }
493 : }
494 : }
495 0 : }
496 :
497 : template <typename TParent>
498 0 : void prolongate(TParent* parent,
499 : TransferValueAccessor& vValueChild,
500 : TransferValueAccessor& vValueParent)
501 : {
502 : const MultiGrid* mg = IElemProlongation<TDomain>::m_spGrid.get();
503 : // a) prolongation from a side
504 0 : if(TParent::dim == m_lfeid.dim()-1){
505 : typedef typename TParent::sideof TElem;
506 : typedef TParent TSide;
507 :
508 : // get child sides
509 0 : const int numChild = mg->num_children<TParent>(parent);
510 0 : if(numChild == 0) return;
511 :
512 : // get childs
513 0 : std::vector<TSide*> vChildSide(numChild);
514 0 : for(int c = 0; c < numChild; ++c)
515 0 : vChildSide[c] = mg->get_child<TParent>(parent, c);
516 :
517 : // get associated macro elements, this elem is a side of
518 : typename Grid::traits<TElem>::secure_container vElem;
519 0 : const_cast<MultiGrid*>(mg)->associated_elements(vElem, parent);
520 0 : std::vector<TElem*> vParentElem(vElem.size());
521 0 : for(size_t p = 0; p < vElem.size(); ++p)
522 0 : vParentElem[p] = vElem[p];
523 :
524 : // call prolongation
525 0 : this->template prolongate<TSide>(vParentElem, vChildSide, vValueChild, vValueParent);
526 0 : }
527 :
528 : // b) prolongation from a element
529 0 : if(TParent::dim == m_lfeid.dim()){
530 : typedef TParent TElem;
531 : typedef typename TParent::side TSide;
532 :
533 : // get child sides
534 0 : const int numChild = mg->num_children<TSide>(parent);
535 0 : if(numChild == 0) return;
536 :
537 : // get childs
538 0 : std::vector<TSide*> vChildSide(numChild);
539 0 : for(int c = 0; c < numChild; ++c)
540 0 : vChildSide[c] = mg->get_child<TSide>(parent, c);
541 :
542 : // get associated macro elements, this elem is a side of
543 : std::vector<TElem*> vParentElem(1);
544 0 : vParentElem[0] = parent;
545 :
546 : // call prolongation
547 0 : this->template prolongate<TSide>(vParentElem, vChildSide, vValueChild, vValueParent);
548 0 : }
549 : };
550 :
551 0 : virtual bool perform_restriction_on(GridBaseObjectId gbo) {
552 : // restriction only on sides
553 0 : if(m_lfeid.dim()-1 == gbo) return true;
554 : return false;
555 : }
556 :
557 : template <typename TParent>
558 0 : void do_restrict(TParent* parent,
559 : TransferValueAccessor& vValueChild,
560 : TransferValueAccessor& vValueParent)
561 : {
562 : const MultiGrid* mg = IElemRestriction<TDomain>::m_spGrid.get();
563 :
564 0 : if(TParent::dim != m_lfeid.dim()-1)
565 0 : UG_THROW("Only restrict on sides.")
566 :
567 : // get child sides
568 0 : const int numChild = mg->num_children<TParent>(parent);
569 0 : if(numChild == 0) return;
570 :
571 : // access the values
572 0 : vValueParent.access_inner(parent);
573 : UG_ASSERT(vValueParent.size() == 1, "Num DoFs per Side must be 1 for CR");
574 :
575 : // reset values
576 0 : vValueParent[0] = 0.0;
577 :
578 : // get childs
579 0 : for(int c = 0; c < numChild; ++c)
580 : {
581 : // get child side
582 0 : TParent* childSide = mg->get_child<TParent>(parent, c);
583 :
584 : // access the value
585 0 : vValueChild.access_inner(childSide);
586 : UG_ASSERT(vValueChild.size() == 1, "Num DoFs per Side must be 1 for CR");
587 :
588 : // sum up values
589 0 : vValueParent[0] += vValueChild[0];
590 : }
591 :
592 0 : vValueParent[0] *= (1./numChild);
593 : }
594 :
595 : protected:
596 : LFEID m_lfeid;
597 : };
598 :
599 :
600 : } // end namespace ug
601 :
602 : #endif /* __H__UG__LIB_DISC__FUNCTION_SPACS__LOCAL_TRANSFER__ */
|