Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_local_ass.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2014: G-CSC, Goethe University Frankfurt
3 * Author: Dmitry Logashenko
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/*
34 * nedelec_local_ass.h
35 * Declarations of the functions for assembling the local stiffness, mass and
36 * weak divergence matrices of the Nedelec-type-1 (Whitney-1) based discretization
37 * of the Maxwell equations, as well as the shape functions and their curls
38 * of the Nedelec element.
39 */
40
41#ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__ROT_ROT_ASS__
42#define __H__UG__PLUGINS__ELECTROMAGNETISM__ROT_ROT_ASS__
43
44// other ug4 modules
45#include "common/common.h"
46#include "lib_grid/lg_base.h"
50
51namespace ug{
52namespace Electromagnetism{
53
56
58
83template <typename TDomain, typename TElem>
85{
86public:
87
89 static const int WDim = TDomain::dim;
90
92 typedef typename TDomain::grid_type grid_type;
93
95 typedef typename TDomain::position_type position_type;
96
99
102
104 static const int dim = ref_elem_type::dim;
105
107 static const size_t numCorners = ref_elem_type::numCorners;
108
110 static const size_t numEdges = ref_elem_type::numEdges;
111
113 static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<WDim>::DimElemList>::maxEdges;
114
115public:
116
119 (
120 const TDomain * domain,
121 TElem * elem,
122 const position_type * corners,
125 )
126 {
127 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
128 };
129
131 static void local_mass
132 (
133 const TDomain * domain,
134 TElem * elem,
135 const position_type * corners,
137 )
138 {
139 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
140 };
141
144 (
145 const TDomain * domain,
146 TElem * elem,
147 const position_type * corners,
148 number B [][numEdges]
149 )
150 {
151 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
152 };
153
155 static void get_shapes
156 (
157 const TDomain * domain,
158 TElem * elem,
159 const position_type * corners,
160 const MathVector<dim> local,
161 MathVector<WDim> shapes []
162 )
163 {
164 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
165 };
166
168 static void interpolate
169 (
170 const TDomain * domain,
171 TElem * elem,
172 const position_type * corners,
173 const number dofs [],
174 const MathVector<dim> local [],
175 const size_t n_pnt,
176 MathVector<WDim> values []
177 )
178 {
179 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
180 };
181
183 static void curl
184 (
185 const TDomain * domain,
186 TElem * elem,
187 const position_type * corners,
188 const number dofs [],
189 MathVector<WDim> & curl_vec
190 )
191 {
192 UG_THROW ("Whitney-1 (Nedelec) shapes not implemented for roid " << ref_elem_type::REFERENCE_OBJECT_ID << ".");
193 };
194};
195
197
254template <typename TDomain, typename TElem>
256{
257public:
258
260 static const int WDim = TDomain::dim;
261
263 typedef typename TDomain::grid_type grid_type;
264
266 typedef typename TDomain::position_type position_type;
267
270
273
275 static const int dim = ref_elem_type::dim;
276
278 static const size_t numCorners = ref_elem_type::numCorners;
279
281 static const size_t numEdges = ref_elem_type::numEdges;
282
284 static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<WDim>::DimElemList>::maxEdges;
285
286private:
287
289 static void compute_W0_grads
290 (
291 const position_type * corners,
293 );
294
296 static void get_edge_corners
297 (
298 const TDomain * domain,
299 TElem * elem,
300 size_t edge_corner [numEdges] [2]
301 );
302
303public:
304
306 static void local_stiffness_and_mass
307 (
308 const TDomain * domain,
309 TElem * elem,
310 const position_type * corners,
313 );
314
316 static void local_mass
317 (
318 const TDomain * domain,
319 TElem * elem,
320 const position_type * corners,
322 );
323
325 static void local_div_matrix
326 (
327 const TDomain * domain,
328 TElem * elem,
329 const position_type * corners,
330 number B [][numEdges]
331 );
332
334 static void get_shapes
335 (
336 const TDomain * domain,
337 TElem * elem,
338 const position_type * corners,
339 const MathVector<dim> local,
340 MathVector<WDim> shapes []
341 );
342
344 static void interpolate
345 (
346 const TDomain * domain,
347 TElem * elem,
348 const position_type * corners,
349 const number dofs [],
350 const MathVector<dim> local [],
351 const size_t n_pnt,
352 MathVector<WDim> values []
353 );
354
356 static void curl
357 (
358 const TDomain * domain,
359 TElem * elem,
360 const position_type * corners,
361 const number dofs [],
362 MathVector<WDim> & curl_vec
363 );
364};
365
367
371template <typename TDomain>
372class NedelecT1_LDisc<TDomain, Triangle> : public NedelecT1_LDisc_forSimplex<TDomain, Triangle>
373{
374};
375
377
381template <typename TDomain>
382class NedelecT1_LDisc<TDomain, Tetrahedron> : public NedelecT1_LDisc_forSimplex<TDomain, Tetrahedron>
383{
384};
385
387
389
397template <typename TDomain, int refDim, int WDim = TDomain::dim>
399{
400public:
401
403 typedef typename TDomain::position_type position_type;
404
406 static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<WDim>::DimElemList>::maxEdges;
407
408public:
410 static void value
411 (
412 const TDomain * domain,
413 GridObject * elem,
414 const position_type corners [],
415 const number dofs [],
416 const MathVector<refDim> local [],
417 const size_t n_pnt,
418 MathVector<WDim> values []
419 )
420 { // This is a generic version: It prints the error message. Cf. the specializations below.
421 UG_THROW ("NedelecInterpolation::value: No implementation of the Nedelec elements for "
422 "Reference Object " << elem->reference_object_id () <<
423 " of reference dim. " << refDim << " in a " << WDim << "d domain.");
424 }
425
427
433 static void curl
434 (
435 const TDomain * domain,
436 GridObject * elem,
437 const position_type corners [],
438 const number dofs [],
439 MathVector<WDim> & curl_vec
440 )
441 { // This is a generic version: It prints the error message. Cf. the specializations below.
442 UG_THROW ("NedelecInterpolation::curl: No implementation of the Nedelec elements for "
443 "Reference Object " << elem->reference_object_id () <<
444 " of reference dim. " << refDim << " in a " << WDim << "d domain.");
445 }
446
448
457 (
458 const TDomain * domain,
459 GridObject * elem,
460 const position_type corners [],
461 const number dofs [],
462 const MathVector<WDim> & normal,
463 const position_type pnt,
464 number & flux
465 )
466 { // This is a generic version: It prints the error message. Cf. the specializations below.
467 UG_THROW ("NedelecInterpolation::curl_flux: No implementation of the Nedelec elements for "
468 "Reference Object " << elem->reference_object_id () <<
469 " of reference dim. " << refDim << " in a " << WDim << "d domain.");
470 return false;
471 }
472};
473
475template <typename TDomain>
476class NedelecInterpolation<TDomain, 2, 2>
477{
478public:
479
481 typedef typename TDomain::position_type position_type;
482
484 static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<2>::DimElemList>::maxEdges;
485
486public:
488 static void value
489 (
490 const TDomain * domain,
491 GridObject * elem,
492 const position_type corners [],
493 const number dofs [],
494 const MathVector<2> local [],
495 const size_t n_pnt,
496 MathVector<2> values []
497 )
498 {
499 if (elem->reference_object_id () != ROID_TRIANGLE)
500 UG_THROW ("No implementation of the Nedelec elements for "
501 "Reference Object " << elem->reference_object_id () <<
502 " of reference dim. 2 in a 2d domain. (This must be a triangle.)");
504 (domain, (Triangle *) elem, corners, dofs, local, n_pnt, values);
505 }
506
508
514 static void curl
515 (
516 const TDomain * domain,
517 GridObject * elem,
518 const position_type corners [],
519 const number dofs [],
520 MathVector<2> & curl_vec
521 )
522 {
523 if (elem->reference_object_id () != ROID_TRIANGLE)
524 UG_THROW ("No implementation of the Nedelec elements for "
525 "Reference Object " << elem->reference_object_id () <<
526 " of reference dim. 2 in a 2d domain. (This must be a triangle.)");
528 (domain, (Triangle *) elem, corners, dofs, curl_vec);
529 }
530
532
541 (
542 const TDomain * domain,
543 GridObject * elem,
544 const position_type corners [],
545 const number dofs [],
546 const MathVector<2> & normal,
547 const position_type pnt,
548 number & flux
549 )
550 { // This is a generic version: It prints the error message. Cf. the specializations below.
551 UG_THROW ("NedelecInterpolation::curl_flux: No implementation of the Nedelec elements for "
552 "Reference Object " << elem->reference_object_id () <<
553 " of reference dim. 3 in a 3d domain.");
554 return false;
555 }
556};
557
559template <typename TDomain>
560class NedelecInterpolation<TDomain, 3, 3>
561{
562public:
563
565 typedef typename TDomain::position_type position_type;
566
568 static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<3>::DimElemList>::maxEdges;
569
570public:
572 static void value
573 (
574 const TDomain * domain,
575 GridObject * elem,
576 const position_type corners [],
577 const number dofs [],
578 const MathVector<3> local [],
579 const size_t n_pnt,
580 MathVector<3> values []
581 )
582 {
584 UG_THROW ("No implementation of the Nedelec elements for "
585 "Reference Object " << elem->reference_object_id () <<
586 " of reference dim. 3 in a 3d domain. (This must be a tetrahedron.)");
588 (domain, (Tetrahedron *) elem, corners, dofs, local, n_pnt, values);
589 }
590
592
597 static void curl
598 (
599 const TDomain * domain,
600 GridObject * elem,
601 const position_type corners [],
602 const number dofs [],
603 MathVector<3> & curl_vec
604 )
605 {
607 UG_THROW ("No implementation of the Nedelec elements for "
608 "Reference Object " << elem->reference_object_id () <<
609 " of reference dim. 3 in a 3d domain. (This must be a tetrahedron.)");
611 (domain, (Tetrahedron *) elem, corners, dofs, curl_vec);
612 }
613
615
623 static number curl_flux
624 (
625 const TDomain * domain,
626 GridObject * elem,
627 const position_type corners [],
628 const number dofs [],
629 const MathVector<3> & normal,
630 const position_type pnt,
631 number & flux
632 );
633};
634
638template <int dim>
639inline void get_ave_vector
640(
641 const size_t n_vec,
642 const MathVector<dim> * vec,
643 MathVector<dim> & ave
644)
645{
646 ave = 0;
647 for (size_t i = 0; i < n_vec; i++) ave += vec [i];
648 ave /= n_vec;
649}
650
651} // end namespace Electromagnetism
652} // end namespace ug
653
655
656#endif /* __H__UG__PLUGINS__ELECTROMAGNETISM__ROT_ROT_ASS__ */
657
658/* End of File */
static void curl(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], MathVector< 2 > &curl_vec)
computes curl of the function
Definition nedelec_local_ass.h:515
static void value(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< 2 > local[], const size_t n_pnt, MathVector< 2 > values[])
computes the values at given points
Definition nedelec_local_ass.h:489
static number curl_flux(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< 2 > &normal, const position_type pnt, number &flux)
computes flux of the curl through a given plane in an element
Definition nedelec_local_ass.h:541
TDomain::position_type position_type
type of the geometric positions (WDim-vectors)
Definition nedelec_local_ass.h:481
TDomain::position_type position_type
type of the geometric positions (WDim-vectors)
Definition nedelec_local_ass.h:565
static void value(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< 3 > local[], const size_t n_pnt, MathVector< 3 > values[])
computes the values at given points
Definition nedelec_local_ass.h:573
static void curl(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], MathVector< 3 > &curl_vec)
computes curl of the function
Definition nedelec_local_ass.h:598
Interpolation of the Nedelec dofs and their curls.
Definition nedelec_local_ass.h:399
TDomain::position_type position_type
type of the geometric positions (WDim-vectors)
Definition nedelec_local_ass.h:403
static number curl_flux(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< WDim > &normal, const position_type pnt, number &flux)
computes flux of the curl through a given plane in an element
Definition nedelec_local_ass.h:457
static void curl(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], MathVector< WDim > &curl_vec)
computes curl of the function
Definition nedelec_local_ass.h:434
static const size_t maxNumEdges
max. number of the edges of the full-dimensional elements in the domain
Definition nedelec_local_ass.h:406
static void value(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< refDim > local[], const size_t n_pnt, MathVector< WDim > values[])
computes the values at given points
Definition nedelec_local_ass.h:411
Helper class for the specialization of NedelecT1_LDisc for simplices (triangles and tetrahedrons)
Definition nedelec_local_ass.h:256
static const size_t maxNumEdges
max. number of the edges of the full-dimensional elements in the domain
Definition nedelec_local_ass.h:284
LagrangeP1< ref_elem_type > W0_shapes_type
shapes and derivatives of the Whitney-0 (Lagrange P1) shape functions
Definition nedelec_local_ass.h:272
static void local_mass(const TDomain *domain, TElem *elem, const position_type *corners, number M[maxNumEdges][maxNumEdges])
assembles the local mass matrix of the Nedelec elements
Definition nedelec_local_ass_impl.h:238
static void interpolate(const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], const MathVector< dim > local[], const size_t n_pnt, MathVector< WDim > values[])
computes of the values of the grid functions
Definition nedelec_local_ass_impl.h:398
TDomain::position_type position_type
type of the geometric positions (WDim-vectors)
Definition nedelec_local_ass.h:266
static const size_t numCorners
total number of the corners
Definition nedelec_local_ass.h:278
static void local_div_matrix(const TDomain *domain, TElem *elem, const position_type *corners, number B[][numEdges])
assembles the discrete weak div operator
Definition nedelec_local_ass_impl.h:321
static void compute_W0_grads(const position_type *corners, MathVector< WDim > grad[numCorners])
computes the gradients of the Whitney-0 (Lagrange P1) shape functions
Definition nedelec_local_ass_impl.h:56
static void get_shapes(const TDomain *domain, TElem *elem, const position_type *corners, const MathVector< dim > local, MathVector< WDim > shapes[])
computes the Nedelec shapes at a given point
Definition nedelec_local_ass_impl.h:362
static void curl(const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], MathVector< WDim > &curl_vec)
computes the curl of the grid functions (in 2d represented as instead of )
Definition nedelec_local_ass_impl.h:449
static const size_t numEdges
total number of the edges
Definition nedelec_local_ass.h:281
static void local_stiffness_and_mass(const TDomain *domain, TElem *elem, const position_type *corners, number S[maxNumEdges][maxNumEdges], number M[maxNumEdges][maxNumEdges])
assembles the local stiffness and mass matrices of the rot-rot operator
Definition nedelec_local_ass_impl.h:145
static void get_edge_corners(const TDomain *domain, TElem *elem, size_t edge_corner[numEdges][2])
gets the correspondence between the edge dof indices and the corners of the element
Definition nedelec_local_ass_impl.h:97
reference_element_traits< TElem >::reference_element_type ref_elem_type
type of reference element
Definition nedelec_local_ass.h:269
TDomain::grid_type grid_type
type of the grid
Definition nedelec_local_ass.h:263
static const int dim
dimension of reference element
Definition nedelec_local_ass.h:275
static const int WDim
world dimention
Definition nedelec_local_ass.h:260
Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators.
Definition nedelec_local_ass.h:85
static void get_shapes(const TDomain *domain, TElem *elem, const position_type *corners, const MathVector< dim > local, MathVector< WDim > shapes[])
computes the Nedelec shapes at a given point
Definition nedelec_local_ass.h:156
static void interpolate(const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], const MathVector< dim > local[], const size_t n_pnt, MathVector< WDim > values[])
computes of the values of the grid functions
Definition nedelec_local_ass.h:169
static void local_stiffness_and_mass(const TDomain *domain, TElem *elem, const position_type *corners, number S[maxNumEdges][maxNumEdges], number M[maxNumEdges][maxNumEdges])
assembles the local stiffness and the mass matrices of the rot-rot operator
Definition nedelec_local_ass.h:119
reference_element_traits< TElem >::reference_element_type ref_elem_type
type of reference element
Definition nedelec_local_ass.h:98
TDomain::grid_type grid_type
type of the grid
Definition nedelec_local_ass.h:92
LagrangeP1< ref_elem_type > W0_shapes_type
shapes and derivatives of the Whitney-0 (Lagrange P1) shape functions
Definition nedelec_local_ass.h:101
static const int dim
dimension of reference element
Definition nedelec_local_ass.h:104
static void local_mass(const TDomain *domain, TElem *elem, const position_type *corners, number M[maxNumEdges][maxNumEdges])
assembles the local mass matrix of the Nedelec elements
Definition nedelec_local_ass.h:132
TDomain::position_type position_type
type of the geometric positions (WDim-vectors)
Definition nedelec_local_ass.h:95
static const size_t maxNumEdges
max. number of the edges of the full-dimensional elements in the domain
Definition nedelec_local_ass.h:113
static void curl(const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], MathVector< WDim > &curl_vec)
computes the curl of the grid functions (in 2d represented as instead of )
Definition nedelec_local_ass.h:184
static const size_t numCorners
total number of the corners
Definition nedelec_local_ass.h:107
static const int WDim
world dimention
Definition nedelec_local_ass.h:89
static void local_div_matrix(const TDomain *domain, TElem *elem, const position_type *corners, number B[][numEdges])
assembles the discrete weak div operator
Definition nedelec_local_ass.h:144
static const size_t numEdges
total number of the edges
Definition nedelec_local_ass.h:110
virtual ReferenceObjectID reference_object_id() const=0
#define UG_THROW(msg)
double number
void get_ave_vector(const size_t n_vec, const MathVector< dim > *vec, MathVector< dim > &ave)
Definition nedelec_local_ass.h:640
ROID_TETRAHEDRON
ROID_TRIANGLE