Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_project.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_project.h - classes for the projection of the functions based
35 * on the Nedelec element to the space of the divergence-free functions.
36 */
37#ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_PROJECT__
38#define __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_PROJECT__
39
40#include "common/common.h"
53
54#ifdef UG_PARALLEL
57#endif
58
59#include "em_material.h"
60
61namespace ug{
62namespace Electromagnetism{
63
65
70template <typename TDomain, typename TAlgebra>
71class NedelecProject : public IPProcessVector<typename TAlgebra::vector_type>
72{
75
76public:
78 typedef TDomain domain_type;
79
81 typedef typename TDomain::grid_type grid_type;
82
84 typedef typename domain_type::subset_handler_type subset_handler_type;
85
87 typedef TAlgebra algebra_type;
89 typedef typename TAlgebra::vector_type vector_type;
91 typedef typename TAlgebra::matrix_type matrix_type;
92
101
103 static const int WDim = TDomain::dim;
104
106 typedef typename TDomain::position_type position_type;
107
108private:
109
111 typedef AInt a_vert_cond_type; // the insulators are marked with -2, not with -1!
113
114public:
115
118 (
119 SmartPtr<EMaterial<TDomain> > emMatherial,
120 SmartPtr<ApproximationSpace<TDomain> > vertApproxSpace,
122 );
123
126 {
127 for (size_t i = 0; i < m_DVF_phi.size (); i++) delete m_DVF_phi [i];
128 }
129
132 (
134 )
135 {
136 m_spDirichlet = spDirichlet;
137 }
138
140 void set_dampDVFs (bool val)
141 {
142 m_bDampDVFs = val;
143 }
144
146 void apply
147 (
148 vector_type & vec
149 );
150
152 void apply
153 (
155 const char * fct_names
156 );
157
159 void apply
160 (
162 const FunctionGroup & fctGrp
163 );
164
166 void compute_div
167 (
169 const char * u_fct_name,
171 );
172
173private:
174
176 inline void project_func
177 (
178 const SmartPtr<TDomain> & domain,
179 const SmartPtr<DoFDistribution> & edgeDD,
180 vector_type & u,
181 size_t fct,
182 const SmartPtr<DoFDistribution> & vertDD,
183 pot_gf_type & aux_rhs,
184 pot_gf_type & aux_cor
185 );
186
187// The weak divergence and gradient operators:
188
190 template <typename TElem>
192 (
193 const TDomain & domain,
194 const DoFDistribution & edgeDD,
195 const vector_type & u,
196 size_t fct,
197 const DoFDistribution & vertDD,
198 pot_vector_type & div
199 );
200
202 struct WeakDiv
203 {
205 const TDomain & m_domain;
208 size_t m_fct;
211
213 (
214 this_type * pThis,
215 const TDomain & domain,
216 const DoFDistribution & edgeDD,
217 const vector_type & u,
218 size_t fct,
219 const DoFDistribution & vertDD,
220 pot_vector_type & div
221 )
222 : m_pThis (pThis), m_domain (domain), m_edgeDD (edgeDD), m_u (u),
223 m_fct (fct), m_vertDD (vertDD), m_div (div) {}
224
225 template <typename TElem> void operator() (TElem &)
226 {
227 m_pThis->template weak_div_elem_type<TElem>
229 }
230 };
231
233 template <typename TElem>
235 (
236 const DoFDistribution & vertDD,
237 pot_vector_type & div,
239 );
240
243 {
248
250 (
251 this_type * pThis,
252 const DoFDistribution & vertDD,
253 pot_vector_type & div,
255 )
256 : m_pThis (pThis), m_vertDD (vertDD), m_div (div), m_charge (charge) {}
257
258 template <typename TElem> void operator() (TElem &)
259 {
260 m_pThis->template clear_div_in_conductors<TElem> (m_vertDD, m_div, m_charge);
261 }
262 };
263
265 void assemble_div
266 (
267 const TDomain & domain,
268 const DoFDistribution & edgeDD,
269 const vector_type & u,
270 size_t fct,
271 const DoFDistribution & vertDD,
272 pot_vector_type & div,
274 );
275
277 void distribute_cor
278 (
279 const TDomain & domain,
280 const DoFDistribution & edgeDD,
281 vector_type & u,
282 size_t fct,
283 const DoFDistribution & vertDD,
284 const pot_vector_type & cor
285 );
286
287// Computation of charges of the DVFs and elimination of the DVFs
288
290 void alloc_DVFs
291 (
292 const TDomain & domain,
293 pot_gf_type & aux_rhs
294 );
295
297 void compute_DVFs
298 (
299 pot_gf_type & aux_rhs
300 );
301
304 (
305 const TDomain & domain,
306 DoFDistribution & vertDD
307 );
308
310 void damp_DVFs
311 (
312 pot_vector_type & cor,
313 const DenseVector<VariableArray1<number> > & charge
314 );
315
317 template <typename TElem>
319 (
320 DoFDistribution & vertDD,
321 aa_vert_cond_type & vert_base_cond
322 );
323
326 {
330
332 (
333 this_type * pThis,
334 DoFDistribution & vertDD,
335 aa_vert_cond_type & vert_base_cond
336 )
337 : m_pThis (pThis), m_vertDD (vertDD), m_vert_base_cond (vert_base_cond) {}
338
339 template <typename TElem> void operator() (TElem &)
340 {
341 m_pThis->template mark_cond_vert_elem_type<TElem> (m_vertDD, m_vert_base_cond);
342 }
343 };
344
346 template <typename TElem>
348 (
349 const TDomain & domain,
350 const DoFDistribution & vertDD,
351 const aa_vert_cond_type & vert_base_cond,
353 );
354
357 {
359 const TDomain & m_domain;
363
365 (
366 this_type * pThis,
367 const TDomain & domain,
368 const DoFDistribution & vertDD,
369 const aa_vert_cond_type & vert_base_cond,
371 )
372 : m_pThis (pThis), m_domain (domain), m_vertDD (vertDD),
373 m_vert_base_cond (vert_base_cond), m_C (C) {}
374
375 template <typename TElem> void operator() (TElem &)
376 {
377 m_pThis->template integrate_div_DVF_elem_type<TElem>
379 }
380 };
381
382// Laplace operators for the auxiliary problems
383
385 template <typename TElem>
387 {
388 public:
390 static const int numCorners = ref_elem_type::numCorners;
391
393 static void stiffness
394 (
395 GridObject * elem,
396 const position_type vCornerCoords [],
397 number loc_A [numCorners] [numCorners]
398 );
399 };
400
402 class AuxLaplaceLocAss : public IElemDisc<TDomain>
403 {
404 private:
406 static const size_t _C_ = 0;
407
410
411 public:
414 (
415 NedelecProject & master
416 );
417
419 virtual void prepare_setting
420 (
421 const std::vector<LFEID> & vLfeID,
422 bool bNonRegular
423 );
424
425 private:
428 template <typename TElem>
429 void prepare_element_loop (ReferenceObjectID roid, int si);
430
431 template <typename TElem>
432 void ass_JA_elem (LocalMatrix & J, const LocalVector & u, GridObject * elem, const position_type vCornerCoords []);
433
434 template <typename TElem>
435 void prepare_element (const LocalVector & u, GridObject * elem, ReferenceObjectID roid, const position_type vCornerCoords [])
436 {};
437 template <typename TElem>
439 {};
440 template <typename TElem>
441 void ass_rhs_elem (LocalVector & d, GridObject * elem, const position_type vCornerCoords [])
442 {};
443 template <typename TElem>
444 void ass_JM_elem (LocalMatrix & J, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
445 {
446 UG_THROW ("NedelecProject: Attempt to use nonstationary discretization for the potential.");
447 };
448 template <typename TElem>
449 void ass_dA_elem (LocalVector & d, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
450 {
451 UG_THROW ("NedelecProject: Attempt to assemble non-linear system for the potential.");
452 };
453 template <typename TElem>
454 void ass_dM_elem (LocalVector & d, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
455 {
456 UG_THROW ("NedelecProject: Attempt to use nonstationary discretization for the potential.");
457 };
459
463
467 template< typename TElem > void operator () (TElem&)
468 {m_pThis->register_loc_discr_func<TElem> ();}
469 };
470
471 template <typename TElem>
474
475 private:
476
479
482 };
483
485 class AuxLaplaceRHS : public IDomainConstraint<TDomain, TPotAlgebra>
486 {
487 private:
490
493
496
497 public:
500 (
501 NedelecProject & master
502 );
503
506 (
507 int base_cond = -1
508 )
509 {
510 m_base_cond = base_cond;
511 }
512
514 void adjust_jacobian
515 (
516 pot_matrix_type & J,
517 const pot_vector_type & u,
519 int type,
520 number time = 0.0,
522 const number s_a0 = 1.0
523 );
524
526 void adjust_defect
527 (
528 pot_vector_type & d,
529 const pot_vector_type & u,
531 int type,
532 number time = 0.0,
534 const std::vector<number> * vScaleMass = NULL,
535 const std::vector<number> * vScaleStiff = NULL
536 );
537
539 void adjust_solution
540 (
541 pot_vector_type & u,
543 int type,
544 number time = 0.0
545 );
546
549 (
550 pot_matrix_type & A,
551 pot_vector_type & b,
553 int type,
554 number time = 0.0
555 )
556 { // Note that this function is not really used, so it needs not to be optimal.
557 AuxLaplaceRHS::adjust_jacobian (A, b, dd, type, time); // the 2nd arg. is dummy: A does not depend on u
559 }
560
563 (
564 pot_vector_type & b,
565 const pot_vector_type & u,
567 int type,
568 number time = 0.0
569 )
570 {
572 }
573
575 int type () const {return CT_DIRICHLET;}
576
577 private:
578
581 (
582 pot_vector_type & u,
583 const DoFDistribution * dd
584 );
585
588 (
589 pot_matrix_type & A,
590 const DoFDistribution * dd
591 );
592
594 template <typename TElem>
596 (
597 int si,
598 number val,
599 pot_vector_type & u,
600 const DoFDistribution * dd
601 );
602
605 {
607 int m_si;
611
613 (
614 AuxLaplaceRHS * pThis,
615 int si,
616 number val,
617 pot_vector_type & u,
618 const DoFDistribution * dd
619 )
620 : m_pThis (pThis), m_si (si), m_val (val), m_u (u), m_dd (dd) {};
621
622 template <typename TElem> void operator() (TElem &)
623 {
624 m_pThis->template set_value_on_subset<TElem> (m_si, m_val, m_u, m_dd);
625 }
626 };
627
629 template <typename TElem>
631 (
632 int si,
633 pot_matrix_type & A,
634 const DoFDistribution * dd
635 );
636
639 {
641 int m_si;
644
646 (
647 AuxLaplaceRHS * pThis,
648 int si,
649 pot_matrix_type & A,
650 const DoFDistribution * dd
651 )
652 : m_pThis (pThis), m_si (si), m_A (A), m_dd (dd) {};
653
654 template <typename TElem> void operator() (TElem &)
655 {
656 m_pThis->template set_identity_on_subset<TElem> (m_si, m_A, m_dd);
657 }
658 };
659
660 private:
661
664
667 };
668
669private:
670
673
676
679
682
694 std::vector<pot_gf_type *> m_DVF_phi;
697};
698
699} // end namespace Electromagnetism
700} // end namespace ug
701
702#include "nedelec_project_impl.h"
703
704#endif // __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_PROJECT__
705
706/* End of File */
Common interface to get the Dirichlet boundary conditions.
Definition em_material.h:287
Class for subdomain-dependent data for the E-based formulated problems.
Definition em_material.h:63
Class for local assembling of the auxiliary Laplace operator.
Definition nedelec_project.h:403
void prepare_element_loop(ReferenceObjectID roid, int si)
Definition nedelec_project_impl.h:926
bool m_do_assemble_here
flag whether to assemble on a given element
Definition nedelec_project.h:481
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegular)
switches between non-regular and regular grids
Definition nedelec_project_impl.h:875
NedelecProject & m_master
the master object of the projection
Definition nedelec_project.h:478
void ass_dM_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_project.h:454
void ass_dA_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_project.h:449
void prepare_element(const LocalVector &u, GridObject *elem, ReferenceObjectID roid, const position_type vCornerCoords[])
Definition nedelec_project.h:435
void register_loc_discr_func()
Definition nedelec_project_impl.h:907
void finish_element_loop()
Definition nedelec_project.h:438
void ass_JA_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
transfer the local stiffness matrix to the global discretization
Definition nedelec_project_impl.h:942
static const size_t _C_
a name for the component in the grid functions
Definition nedelec_project.h:406
void register_all_loc_discr_funcs()
Definition nedelec_project_impl.h:895
void ass_JM_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_project.h:444
void ass_rhs_elem(LocalVector &d, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_project.h:441
static const size_t maxCorners
maximum number of corners of an element
Definition nedelec_project.h:409
Constraint that assembles the rhs und the bc for the auxiliary problems.
Definition nedelec_project.h:486
void adjust_defect(pot_vector_type &d, const pot_vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< pot_vector_type > > vSol=SPNULL, const std::vector< number > *vScaleMass=NULL, const std::vector< number > *vScaleStiff=NULL)
sets a zero value in the defect for all conductor indices
Definition nedelec_project_impl.h:1171
void set_base_conductor(int base_cond=-1)
sets the base conductor index (or -1 for the divergence correction)
Definition nedelec_project.h:506
int type() const
returns the type of the constraints
Definition nedelec_project.h:575
void adjust_linear(pot_matrix_type &A, pot_vector_type &b, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets unity rows in A and dirichlet values in right-hand side b
Definition nedelec_project.h:549
void set_value_on_subset(int si, number val, pot_vector_type &u, const DoFDistribution *dd)
sets constant value on the closure of a full-dimensional subset
Definition nedelec_project_impl.h:1070
void set_identity_Dirichlet(pot_matrix_type &A, const DoFDistribution *dd)
sets identity matrix on the whole Dirichlet boundary
Definition nedelec_project_impl.h:1025
NedelecProject & m_master
the master object of the projection
Definition nedelec_project.h:663
void set_identity_on_subset(int si, pot_matrix_type &A, const DoFDistribution *dd)
sets identity matrix on the closure of a full-dimensional subset
Definition nedelec_project_impl.h:1108
DoFDistribution::traits< Vertex >::const_iterator t_vert_iterator
Iterator over vertices.
Definition nedelec_project.h:495
int m_base_cond
the base conductor index (or -1 for insulators)
Definition nedelec_project.h:666
static const size_t maxCorners
maximum number of corners of an element
Definition nedelec_project.h:492
void adjust_rhs(pot_vector_type &b, const pot_vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets the dirichlet value in the right-hand side
Definition nedelec_project.h:563
void set_zero_Dirichlet(pot_vector_type &u, const DoFDistribution *dd)
sets zero values on the whole Dirichlet boundary
Definition nedelec_project_impl.h:981
void adjust_jacobian(pot_matrix_type &J, const pot_vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< pot_vector_type > > vSol=SPNULL, const number s_a0=1.0)
sets a unity row for all conductor indices
Definition nedelec_project_impl.h:1144
DoFDistribution::traits< Edge >::const_iterator t_edge_iterator
Iterator over edges.
Definition nedelec_project.h:489
void adjust_solution(pot_vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets the value in the solution for all conductor indices
Definition nedelec_project_impl.h:1199
Helper class for assembling the local stiffness matrix of Laplacian in various routines.
Definition nedelec_project.h:387
static const int numCorners
Definition nedelec_project.h:390
static void stiffness(GridObject *elem, const position_type vCornerCoords[], number loc_A[numCorners][numCorners])
Computes the local discretization of the Laplace operator.
Definition nedelec_project_impl.h:828
reference_element_traits< TElem >::reference_element_type ref_elem_type
Definition nedelec_project.h:389
Projection procedure class.
Definition nedelec_project.h:72
void project_func(const SmartPtr< TDomain > &domain, const SmartPtr< DoFDistribution > &edgeDD, vector_type &u, size_t fct, const SmartPtr< DoFDistribution > &vertDD, pot_gf_type &aux_rhs, pot_gf_type &aux_cor)
Projects one function (i.e. performs the main action):
Definition nedelec_project_impl.h:247
void mark_cond_vert_elem_type(DoFDistribution &vertDD, aa_vert_cond_type &vert_base_cond)
Sets the base conductor indices for every vertex (for all elements of the same type)
Definition nedelec_project_impl.h:702
void assemble_div(const TDomain &domain, const DoFDistribution &edgeDD, const vector_type &u, size_t fct, const DoFDistribution &vertDD, pot_vector_type &div, DenseVector< VariableArray1< number > > &charge)
Assembles the weak divergence operator on the whole grid.
Definition nedelec_project_impl.h:404
AInt a_vert_cond_type
Type of the attachment and its accessor the base conductor indices.
Definition nedelec_project.h:111
void distribute_cor(const TDomain &domain, const DoFDistribution &edgeDD, vector_type &u, size_t fct, const DoFDistribution &vertDD, const pot_vector_type &cor)
Updates the grid function by the potential correction.
Definition nedelec_project_impl.h:461
TPotAlgebra::vector_type pot_vector_type
Vector type for the potential space.
Definition nedelec_project.h:96
SmartPtr< AuxLaplaceRHS > m_auxLaplaceRHS
Rhs and BC for the auxiliary problem.
Definition nedelec_project.h:686
void set_Dirichlet(SmartPtr< EMDirichlet< TDomain, TAlgebra > > spDirichlet)
Sets the Dirichlet boundary.
Definition nedelec_project.h:132
SmartPtr< AssembledLinearOperator< TPotAlgebra > > m_auxLaplaceOp
Matrix of the discretization.
Definition nedelec_project.h:690
bool m_bDampDVFs
Whether to damp the Dirichlet vector fields (DVFs)
Definition nedelec_project.h:678
Grid::VertexAttachmentAccessor< a_vert_cond_type > aa_vert_cond_type
Definition nedelec_project.h:112
GridFunction< TDomain, TPotAlgebra > pot_gf_type
Grid function type for the potential space.
Definition nedelec_project.h:100
TAlgebra algebra_type
Type of algebra (for the Nedelec-element-based grid functions)
Definition nedelec_project.h:87
TDomain domain_type
Type of Domain.
Definition nedelec_project.h:78
~NedelecProject()
Destructor.
Definition nedelec_project.h:125
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
Approximation space for the potential (vertex-centered) space.
Definition nedelec_project.h:675
TPotAlgebra::matrix_type pot_matrix_type
Matrix type for the potential space.
Definition nedelec_project.h:98
void compute_DVF_potential_coeffs(const TDomain &domain, DoFDistribution &vertDD)
Computes the potential coefficients.
Definition nedelec_project_impl.h:617
TAlgebra::vector_type vector_type
Type of Vector (for the Nedelec-element-based grid functions)
Definition nedelec_project.h:89
SmartPtr< EMDirichlet< TDomain, TAlgebra > > m_spDirichlet
Dirichlet boundary.
Definition nedelec_project.h:681
void set_dampDVFs(bool val)
Sets whether to damp Dirichlet vector fields.
Definition nedelec_project.h:140
void clear_div_in_conductors(const DoFDistribution &vertDD, pot_vector_type &div, DenseVector< VariableArray1< number > > &charge)
Sets the div operator to 0 in conductors (for all elements of the same type)
Definition nedelec_project_impl.h:358
void compute_div(SmartPtr< GridFunction< TDomain, TAlgebra > > sp_u, const char *u_fct_name, SmartPtr< GridFunction< TDomain, TPotAlgebra > > sp_div)
Computes the weak divergence in insulators.
Definition nedelec_project_impl.h:192
SmartPtr< DomainDiscretization< TDomain, TPotAlgebra > > m_auxLaplaceAss
Global discretization of the auxiliary equations.
Definition nedelec_project.h:688
SmartPtr< ILinearOperatorInverse< pot_vector_type > > m_potSolver
Solver for the auxliliary equations.
Definition nedelec_project.h:692
domain_type::subset_handler_type subset_handler_type
Type of subset handler.
Definition nedelec_project.h:84
SmartPtr< AuxLaplaceLocAss > m_auxLocLaplace
Local discretization of the auxiliary equations.
Definition nedelec_project.h:684
void apply(vector_type &vec)
Performs the projection of all functions in a grid function.
Definition nedelec_project_impl.h:90
DenseMatrix< VariableArray2< number > > m_potCoeff
Potential coefficients:
Definition nedelec_project.h:696
NedelecProject< TDomain, TAlgebra > this_type
This type.
Definition nedelec_project.h:74
static const int WDim
world dimention
Definition nedelec_project.h:103
void alloc_DVFs(const TDomain &domain, pot_gf_type &aux_rhs)
Allocates memory for the DVFs associated with the ungrounded conductors.
Definition nedelec_project_impl.h:516
std::vector< pot_gf_type * > m_DVF_phi
Dirichlet vector fields (DVFs)
Definition nedelec_project.h:694
TDomain::position_type position_type
position type
Definition nedelec_project.h:106
void damp_DVFs(pot_vector_type &cor, const DenseVector< VariableArray1< number > > &charge)
Damps the Dirichlet vector fields (DVFs)
Definition nedelec_project_impl.h:680
void integrate_div_DVF_elem_type(const TDomain &domain, const DoFDistribution &vertDD, const aa_vert_cond_type &vert_base_cond, DenseMatrix< VariableArray2< number > > &C)
Integrates div E over boundaries of conductors for elements of the same type.
Definition nedelec_project_impl.h:746
TDomain::grid_type grid_type
Type of Grid:
Definition nedelec_project.h:81
void compute_DVFs(pot_gf_type &aux_rhs)
Computes the Dirichlet vector fields (DVFs)
Definition nedelec_project_impl.h:586
SmartPtr< EMaterial< TDomain > > m_spEmMaterial
Properties of the matherials in the domain.
Definition nedelec_project.h:672
void weak_div_elem_type(const TDomain &domain, const DoFDistribution &edgeDD, const vector_type &u, size_t fct, const DoFDistribution &vertDD, pot_vector_type &div)
Computes the weak divergence (for all elements of the same type)
Definition nedelec_project_impl.h:289
TAlgebra::matrix_type matrix_type
Type of Vector (for the Nedelec-element-based grid functions)
Definition nedelec_project.h:91
CPUAlgebra TPotAlgebra
The auxiliary algebra type for the space of the potential. (Note: It should be scalar....
Definition nedelec_project.h:94
ConstSmartPtr< DoFDistribution > dd(const GridLevel &gl) const
#define UG_THROW(msg)
double number
CT_DIRICHLET
ReferenceObjectID
AuxLaplaceLocAss * m_pThis
Definition nedelec_project.h:466
RegisterLocalDiscr(AuxLaplaceLocAss *pThis)
Definition nedelec_project.h:465
Helper class for 'set_value_on_subset'.
Definition nedelec_project.h:639
SetIdentityOnSubset(AuxLaplaceRHS *pThis, int si, pot_matrix_type &A, const DoFDistribution *dd)
Definition nedelec_project.h:646
AuxLaplaceRHS * m_pThis
Definition nedelec_project.h:640
const DoFDistribution * m_dd
Definition nedelec_project.h:643
pot_matrix_type & m_A
Definition nedelec_project.h:642
void operator()(TElem &)
Definition nedelec_project.h:654
Helper class for 'set_value_on_subset'.
Definition nedelec_project.h:605
SetValueOnSubset(AuxLaplaceRHS *pThis, int si, number val, pot_vector_type &u, const DoFDistribution *dd)
Definition nedelec_project.h:613
void operator()(TElem &)
Definition nedelec_project.h:622
pot_vector_type & m_u
Definition nedelec_project.h:609
const DoFDistribution * m_dd
Definition nedelec_project.h:610
AuxLaplaceRHS * m_pThis
Definition nedelec_project.h:606
Helper class for assembling the weak divergence operator.
Definition nedelec_project.h:243
this_type * m_pThis
Definition nedelec_project.h:244
pot_vector_type & m_div
Definition nedelec_project.h:246
void operator()(TElem &)
Definition nedelec_project.h:258
ClearDivInConductors(this_type *pThis, const DoFDistribution &vertDD, pot_vector_type &div, DenseVector< VariableArray1< number > > &charge)
Definition nedelec_project.h:250
const DoFDistribution & m_vertDD
Definition nedelec_project.h:245
DenseVector< VariableArray1< number > > & m_charge
Definition nedelec_project.h:247
Helper class for computation of the charges of the DVFs.
Definition nedelec_project.h:357
IntegrateDivDVF(this_type *pThis, const TDomain &domain, const DoFDistribution &vertDD, const aa_vert_cond_type &vert_base_cond, DenseMatrix< VariableArray2< number > > &C)
Definition nedelec_project.h:365
const TDomain & m_domain
Definition nedelec_project.h:359
const aa_vert_cond_type & m_vert_base_cond
Definition nedelec_project.h:361
const DoFDistribution & m_vertDD
Definition nedelec_project.h:360
DenseMatrix< VariableArray2< number > > & m_C
Definition nedelec_project.h:362
void operator()(TElem &)
Definition nedelec_project.h:375
this_type * m_pThis
Definition nedelec_project.h:358
Helper class for setting the base conductor indices to vertices.
Definition nedelec_project.h:326
DoFDistribution & m_vertDD
Definition nedelec_project.h:328
void operator()(TElem &)
Definition nedelec_project.h:339
MarkCondVert(this_type *pThis, DoFDistribution &vertDD, aa_vert_cond_type &vert_base_cond)
Definition nedelec_project.h:332
this_type * m_pThis
Definition nedelec_project.h:327
aa_vert_cond_type & m_vert_base_cond
Definition nedelec_project.h:329
Helper class for the computation the weak divergence.
Definition nedelec_project.h:203
pot_vector_type & m_div
Definition nedelec_project.h:210
void operator()(TElem &)
Definition nedelec_project.h:225
const TDomain & m_domain
Definition nedelec_project.h:205
const DoFDistribution & m_edgeDD
Definition nedelec_project.h:206
size_t m_fct
Definition nedelec_project.h:208
this_type * m_pThis
Definition nedelec_project.h:204
WeakDiv(this_type *pThis, const TDomain &domain, const DoFDistribution &edgeDD, const vector_type &u, size_t fct, const DoFDistribution &vertDD, pot_vector_type &div)
Definition nedelec_project.h:213
const DoFDistribution & m_vertDD
Definition nedelec_project.h:209
const vector_type & m_u
Definition nedelec_project.h:207