Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
local_transfer.h
Go to the documentation of this file.
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
41
42namespace ug{
43
44
45template <typename TDomain>
47 : public ElemProlongationBase<TDomain, PiecewiseConstantElemTransfer<TDomain> >
48 , public ElemRestrictionBase<TDomain, PiecewiseConstantElemTransfer<TDomain> >
49{
50 public:
52
54 if(m_lfeid.dim() == gbo) return true;
55 return false;
56 }
57
58 template <typename TElem>
59 void prolongate(TElem* parent,
60 TransferValueAccessor& vValueChild,
61 TransferValueAccessor& vValueParent)
62 {
64 const int numChild = mg->num_children<TElem>(parent);
65
66 vValueParent.access_inner(parent);
67 for(int c = 0; c < numChild; ++c){
68 TElem* child = mg->get_child<TElem>(parent, c);
69
70 vValueChild.access_inner(child);
71 vValueChild[0] = vValueParent[0];
72 }
73 }
74
76 if(m_lfeid.dim() == gbo) return true;
77 return false;
78 }
79
80 template <typename TElem>
81 void do_restrict(TElem* parent,
82 TransferValueAccessor& vValueChild,
83 TransferValueAccessor& vValueParent)
84 {
86 const int numChild = mg->num_children<TElem>(parent);
87
88 vValueParent.access_inner(parent);
89 UG_ASSERT(vValueParent.size()==1, "Exactly one DoF for Piecewise-Constant")
90 vValueParent[0] = 0.0;
91
92 for(int c = 0; c < numChild; ++c){
93 TElem* child = mg->get_child<TElem>(parent, c);
94
95 vValueChild.access_inner(child);
96 UG_ASSERT(vValueChild.size()==1, "Exactly one DoF for Piecewise-Constant")
97 vValueParent[0] += vValueChild[0];
98 }
99 vValueParent[0] *= (1./numChild);
100 }
101
102 protected:
104};
105
106template <typename TDomain>
108 : public ElemProlongationBase<TDomain, P1LagrangeElemTransfer<TDomain> >
109 , public ElemRestrictionBase<TDomain, P1LagrangeElemTransfer<TDomain> >
110{
111 public:
112 P1LagrangeElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
113
115 if(m_lfeid.dim() < gbo) return false;
116 return true;
117 }
118
119 template <typename TElem>
120 void prolongate(TElem* parent,
121 TransferValueAccessor& vValueChild,
122 TransferValueAccessor& vValueParent)
123 {
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 vValueChild.access_inner(child);
132 vValueParent.access_closure(parent);
133
134 if (vValueChild.size() == 0)
135 return; // current fct not defined on child
136
137 const int numVrt = vValueParent.size();
138 UG_ASSERT(numVrt > 0, "No function value found on parent vertices.");
139
140 vValueChild[0] = vValueParent[0];
141 for(int i = 1; i < numVrt; ++i)
142 vValueChild[0] += vValueParent[i];
143 vValueChild[0] *= 1.0/numVrt;
144 }
145
147 if(gbo == VERTEX) return true;
148 return false;
149 }
150
151 // the following line silences -Woverloaded-virtual
154 TransferValueAccessor& vValueChild,
155 TransferValueAccessor& vValueParent)
156 {
157 UG_THROW("Should not be called.")
158 }
159
160 void do_restrict(Vertex* parent,
161 TransferValueAccessor& vValueChild,
162 TransferValueAccessor& vValueParent)
163 {
165 const int numChild = mg->num_children<Vertex>(parent);
166 if(numChild != 1) UG_THROW("Num child Vertex must be 1");
167
168 Vertex* child = mg->get_child<Vertex>(parent, 0);
169
170 vValueChild.access_inner(child);
171 vValueParent.access_closure(parent);
172
173 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 vValueParent[0] = vValueChild[0];
179 }
180
181 protected:
183};
184
185
186template <typename TDomain>
188 : public ElemProlongationBase<TDomain,StdLagrangeElemTransfer<TDomain> >
189 , public ElemRestrictionBase<TDomain,StdLagrangeElemTransfer<TDomain> >
190{
191 public:
193 static const int dim = TDomain::dim;
194
195 public:
196 StdLagrangeElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
197
199 if(m_lfeid.order() == 1 && gbo != VERTEX) return false;
200 if(m_lfeid.dim() < gbo) return false;
201 return true;
202 }
203
204
205 template <typename TParent, typename TChild>
206 void prolongate(TParent* parent,
207 TransferValueAccessor& vValueChild,
208 TransferValueAccessor& vValueParent,
210 const std::vector<MathVector<dim> >& vCornerParent,
211 const ReferenceObjectID parentRoid)
212 {
214 static const int pdim = TParent::dim;
215
216 // number of children of the base type
217 const int numChild = mg->num_children<TChild>(parent);
218
219 // loop all children
220 for(int c = 0; c < numChild; ++c)
221 {
222 // get child
223 TChild* child = mg->get_child<TChild>(parent, c);
224
225 // access the values
226 vValueChild.access_inner(child);
227
228 // global positions of fine dofs
229 std::vector<MathVector<dim> > vDoFPos;
231
232 // get Reference Mapping
234 ReferenceMappingProvider::get<pdim, dim>(parentRoid, vCornerParent);
235
236 // get local position of DoF
237 std::vector<MathVector<pdim> > vLocPos(vDoFPos.size(), 0.0);
238 map.global_to_local(vLocPos, vDoFPos);
239
240 // evaluate coarse shape fct at fine local point
241 std::vector<std::vector<number> > vvShape;
242 lsfs.shapes(vvShape, vLocPos);
243
244 // sum up interpolated values
245 for(size_t csh = 0; csh < vValueChild.size(); ++csh){
246 vValueChild[csh] = 0.0;
247 for(size_t psh = 0; psh < vValueParent.size(); ++psh)
248 vValueChild[csh] += vvShape[csh][psh] * vValueParent[psh];
249 }
250 }
251 }
252
253
254 template <typename TElem>
255 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 const ReferenceObjectID parentRoid = parent->reference_object_id();
263
264 // access the values
265 vValueParent.access_closure(parent);
266
267 // get vertices if required by some prolongation
268 std::vector<MathVector<dim> > vCornerParent;
270
271 // get local finite element trial spaces
272 const LocalShapeFunctionSet<pdim>& lsfs =
273 LocalFiniteElementProvider::get<pdim>(parentRoid, m_lfeid);
274
275 // prolongate to all children with dimension pdim and lower
276 switch(pdim){
277 case 3: prolongate<TElem,Volume>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
278 case 2: prolongate<TElem,Face>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
279 case 1: prolongate<TElem,Edge>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
280 prolongate<TElem,Vertex>(parent, vValueChild, vValueParent, lsfs, vCornerParent, parentRoid);
281 break;
282 default: UG_THROW("Dimension "<<pdim<<" not supported");
283 }
284 }
285
286 virtual void prolongate(Vertex* parent)
287 {
291
292 // access the values
293 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 vValueChild.access_inner(child);
305
306 // copy value
307 vValueChild[0] = vValueParent[0];
308 }
309
310
312 if(m_lfeid.dim() < gbo) return false;
313 return true;
314 }
315
316 template <typename TElem>
317 void do_restrict(TElem* parent,
318 TransferValueAccessor& vValueChild,
319 TransferValueAccessor& vValueParent)
320 {
321 // access the values
322 vValueParent.access_inner(parent);
323
324 // global positions of fine dofs
325 std::vector<MathVector<dim> > vDoFPos;
327
329
330 // number of children of the base type
331 const int numChild = mg->num_children<TElem>(parent);
332
333 // loop all children
334 for(size_t ip = 0; ip < vDoFPos.size(); ++ip)
335 {
336 int c = 0;
337 for(; c < numChild; ++c)
338 {
339 // get child
340 TElem* child = mg->get_child<TElem>(parent, c);
341
342 // global positions of fine dofs
343 std::vector<MathVector<dim> > vChildDoFPos;
345
346 // check if pos matches
347 int cip = -1;
348 for(size_t i = 0; i < vChildDoFPos.size(); ++i)
349 if(VecDistanceSq(vChildDoFPos[i], vDoFPos[ip]) < 1e-8)
350 cip = i;
351
352 // if not found in this child, continue
353 if(cip == -1) continue;
354
355 // access the values
356 vValueChild.access_closure(child);
357
358 // sum up interpolated values
359 vValueParent[ip] = vValueChild[cip];
360
361 break;
362 }
363 // check that DoF found
364 if(c == numChild)
365 UG_THROW("Lagrange-Element-Coarsening: Cannot find DoF position "
366 "in child elems for "<<ip<<"'s Pos: "<<vDoFPos[ip])
367 }
368
369 }
370
371
372 virtual void do_restrict(Vertex* parent)
373 {
375 const int numChild = mg->num_children<Vertex>(parent);
376 if(numChild != 1) UG_THROW("Num child Vertex must be 1");
377
380
381 Vertex* child = mg->get_child<Vertex>(parent, 0);
382
383 vValueChild.access_inner(child);
384 vValueParent.access_inner(parent);
385
386 vValueParent[0] = vValueChild[0];
387 }
388
389
390 protected:
392};
393
394
395template <typename TDomain>
397 : public ElemProlongationBase<TDomain, CrouzeixRaviartElemTransfer<TDomain> >
398 , public ElemRestrictionBase<TDomain, CrouzeixRaviartElemTransfer<TDomain> >
399{
400 public:
402 static const int dim = TDomain::dim;
403
404 public:
405 CrouzeixRaviartElemTransfer(const LFEID& lfeid) : m_lfeid(lfeid) {}
406
408 // prolongation from elems that have a side as child
409 if(m_lfeid.dim()-1 == gbo) return true;
410 if(m_lfeid.dim() == gbo) return true;
411 return false;
412 }
413
414 // the following line silences -Woverloaded-virtual
416 void prolongate(Vertex* parent,
417 TransferValueAccessor& vValueChild,
418 TransferValueAccessor& vValueParent)
419 {
420 UG_THROW("No implementation for Vertex and Crouzeix-Raviart.")
421 }
422
423 template <typename TSide>
424 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 if(vParentElem.size() > 2)
432 UG_THROW("More than 2 Elements share a side.")
433
434 // reset values to zero
435 for(size_t c = 0; c < vChildSide.size(); ++c){
436 vValueChild.access_inner(vChildSide[c]);
437 for(size_t csh = 0; csh < vValueChild.size(); ++csh)
438 vValueChild[csh] = 0.0;
439 }
440
441 // loop parent Elems
442 for(size_t p = 0; p < vParentElem.size(); ++p)
443 {
444 // get parent elem
445 TElem* parentElem = vParentElem[p];
446
447 // get reference object ids
448 static const int pdim = TElem::dim;
449 const ReferenceObjectID parentRoid = parentElem->reference_object_id();
450
451 // access the values
452 vValueParent.access_closure(parentElem);
453
454 // get vertices if required by some prolongation
455 std::vector<MathVector<dim> > vCornerParent;
457
458 // get local finite element trial spaces
459 const LocalShapeFunctionSet<pdim>& lsfs =
460 LocalFiniteElementProvider::get<pdim>(parentRoid, m_lfeid);
461
462 // loop child sides
463 for(size_t c = 0; c < vChildSide.size(); ++c)
464 {
465 // get child side
466 TSide* childSide = vChildSide[c];
467
468 // access the value
469 vValueChild.access_inner(childSide);
470
471 // global positions of fine dofs
472 std::vector<MathVector<dim> > vDoFPos;
474
475 // get Reference Mapping
477 ReferenceMappingProvider::get<pdim, dim>(parentRoid, vCornerParent);
478
479 // get local position of DoF
480 std::vector<MathVector<pdim> > vLocPos(vDoFPos.size(), 0.0);
481 map.global_to_local(vLocPos, vDoFPos);
482
483 // evaluate coarse shape fct at fine local point
484 std::vector<std::vector<number> > vvShape;
485 lsfs.shapes(vvShape, vLocPos);
486
487 // sum up interpolated values
488 for(size_t csh = 0; csh < vValueChild.size(); ++csh){
489 for(size_t psh = 0; psh < vValueParent.size(); ++psh)
490 vValueChild[csh] += (1./vParentElem.size()) *
491 vvShape[csh][psh] * vValueParent[psh];
492 }
493 }
494 }
495 }
496
497 template <typename TParent>
498 void prolongate(TParent* parent,
499 TransferValueAccessor& vValueChild,
500 TransferValueAccessor& vValueParent)
501 {
503 // a) prolongation from a side
504 if(TParent::dim == m_lfeid.dim()-1){
505 typedef typename TParent::sideof TElem;
506 typedef TParent TSide;
507
508 // get child sides
509 const int numChild = mg->num_children<TParent>(parent);
510 if(numChild == 0) return;
511
512 // get childs
513 std::vector<TSide*> vChildSide(numChild);
514 for(int c = 0; c < numChild; ++c)
515 vChildSide[c] = mg->get_child<TParent>(parent, c);
516
517 // get associated macro elements, this elem is a side of
519 const_cast<MultiGrid*>(mg)->associated_elements(vElem, parent);
520 std::vector<TElem*> vParentElem(vElem.size());
521 for(size_t p = 0; p < vElem.size(); ++p)
522 vParentElem[p] = vElem[p];
523
524 // call prolongation
525 this->template prolongate<TSide>(vParentElem, vChildSide, vValueChild, vValueParent);
526 }
527
528 // b) prolongation from a element
529 if(TParent::dim == m_lfeid.dim()){
530 typedef TParent TElem;
531 typedef typename TParent::side TSide;
532
533 // get child sides
534 const int numChild = mg->num_children<TSide>(parent);
535 if(numChild == 0) return;
536
537 // get childs
538 std::vector<TSide*> vChildSide(numChild);
539 for(int c = 0; c < numChild; ++c)
540 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 vParentElem[0] = parent;
545
546 // call prolongation
547 this->template prolongate<TSide>(vParentElem, vChildSide, vValueChild, vValueParent);
548 }
549 };
550
552 // restriction only on sides
553 if(m_lfeid.dim()-1 == gbo) return true;
554 return false;
555 }
556
557 template <typename TParent>
558 void do_restrict(TParent* parent,
559 TransferValueAccessor& vValueChild,
560 TransferValueAccessor& vValueParent)
561 {
563
564 if(TParent::dim != m_lfeid.dim()-1)
565 UG_THROW("Only restrict on sides.")
566
567 // get child sides
568 const int numChild = mg->num_children<TParent>(parent);
569 if(numChild == 0) return;
570
571 // access the values
572 vValueParent.access_inner(parent);
573 UG_ASSERT(vValueParent.size() == 1, "Num DoFs per Side must be 1 for CR");
574
575 // reset values
576 vValueParent[0] = 0.0;
577
578 // get childs
579 for(int c = 0; c < numChild; ++c)
580 {
581 // get child side
582 TParent* childSide = mg->get_child<TParent>(parent, c);
583
584 // access the value
585 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 vValueParent[0] += vValueChild[0];
590 }
591
592 vValueParent[0] *= (1./numChild);
593 }
594
595 protected:
597};
598
599
600} // end namespace ug
601
602#endif /* __H__UG__LIB_DISC__FUNCTION_SPACS__LOCAL_TRANSFER__ */
parameterString p
Definition local_transfer.h:399
CrouzeixRaviartElemTransfer(const LFEID &lfeid)
Definition local_transfer.h:405
void do_restrict(TParent *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:558
virtual bool perform_restriction_on(GridBaseObjectId gbo)
Definition local_transfer.h:551
void prolongate(const std::vector< typename TSide::sideof * > &vParentElem, const std::vector< TSide * > &vChildSide, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:424
virtual bool perform_prolongation_on(GridBaseObjectId gbo)
Definition local_transfer.h:407
void prolongate(Vertex *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:416
LFEID m_lfeid
Definition local_transfer.h:596
static const int dim
world dimension
Definition local_transfer.h:402
void prolongate(TParent *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:498
virtual base class for reference mappings
Definition reference_mapping_provider.h:53
virtual void global_to_local(MathVector< dim > &locPos, const MathVector< worldDim > &globPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate
Definition local_transfer_interface.h:108
Definition local_transfer_interface.h:175
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
Definition local_transfer_interface.h:76
Definition local_transfer_interface.h:143
Identifier for Local Finite Elements.
Definition local_finite_element_id.h:98
int order() const
returns the order of the local finite element
Definition local_finite_element_id.h:126
int dim() const
returns the space dimension of the local finite element
Definition local_finite_element_id.h:129
virtual base class for local shape function sets
Definition local_shape_function_set.h:70
virtual void shapes(shape_type *vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition multi_grid.h:72
TChild * get_child(TElem *elem, size_t ind) const
returns the i-th child of the given child-type
Definition multi_grid.h:268
size_t num_children(TElem *elem) const
returns the number of children of the given child-type
Definition multi_grid.h:225
Definition local_transfer.h:110
P1LagrangeElemTransfer(const LFEID &lfeid)
Definition local_transfer.h:112
void do_restrict(GridObject *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:153
LFEID m_lfeid
Definition local_transfer.h:182
virtual bool perform_restriction_on(GridBaseObjectId gbo)
Definition local_transfer.h:146
virtual bool perform_prolongation_on(GridBaseObjectId gbo)
Definition local_transfer.h:114
void do_restrict(Vertex *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:160
void prolongate(TElem *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:120
Definition local_transfer.h:49
virtual bool perform_restriction_on(GridBaseObjectId gbo)
Definition local_transfer.h:75
LFEID m_lfeid
Definition local_transfer.h:103
void do_restrict(TElem *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:81
void prolongate(TElem *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:59
PiecewiseConstantElemTransfer(const LFEID &lfeid)
Definition local_transfer.h:51
virtual bool perform_prolongation_on(GridBaseObjectId gbo)
Definition local_transfer.h:53
Definition local_transfer.h:190
static const int dim
world dimension
Definition local_transfer.h:193
void prolongate(TElem *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:255
LFEID m_lfeid
Definition local_transfer.h:391
virtual bool perform_restriction_on(GridBaseObjectId gbo)
Definition local_transfer.h:311
virtual bool perform_prolongation_on(GridBaseObjectId gbo)
Definition local_transfer.h:198
virtual void do_restrict(Vertex *parent)
Definition local_transfer.h:372
virtual void prolongate(Vertex *parent)
Definition local_transfer.h:286
StdLagrangeElemTransfer(const LFEID &lfeid)
Definition local_transfer.h:196
void prolongate(TParent *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent, const LocalShapeFunctionSet< TParent::dim > &lsfs, const std::vector< MathVector< dim > > &vCornerParent, const ReferenceObjectID parentRoid)
Definition local_transfer.h:206
void do_restrict(TElem *parent, TransferValueAccessor &vValueChild, TransferValueAccessor &vValueParent)
Definition local_transfer.h:317
Definition local_transfer_interface.h:46
virtual void access_inner(GridObject *elem)=0
size_t size() const
Definition local_transfer_interface.h:64
virtual void access_closure(GridObject *elem)=0
Base-class for all vertex-types.
Definition grid_base_objects.h:231
void CollectCornerCoordinates(std::vector< typename TAAPos::ValueType > &vCornerCoordsOut, const TElem &elem, const TAAPos &aaPos, bool clearContainer=true)
returns the corner coordinates of a geometric object
Definition domain_util_impl.h:75
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_THROW(msg)
Definition error.h:57
vector_t::value_type VecDistanceSq(const vector_t &v1, const vector_t &v2)
returns the squared distance of two vector_ts.
Definition math_vector_functions_common_impl.hpp:351
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
bool DoFPosition(std::vector< MathVector< dim > > &vPos, const ReferenceObjectID roid, const std::vector< MathVector< dim > > &vCornerCoord, const LFEID &lfeID)
Definition dof_position_util.cpp:202
GridBaseObjectId
enumeration of the GridBaseObjects that make up a grid.
Definition grid_base_objects.h:59
@ VERTEX
Definition grid_base_objects.h:60
bool InnerDoFPosition(std::vector< MathVector< dim > > &vPos, const ReferenceObjectID roid, const std::vector< MathVector< dim > > &vCornerCoord, const LFEID &lfeID)
Definition dof_position_util.cpp:112
PointerConstArray< TElem * > secure_container
Definition grid.h:146