Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_source.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
37#ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_SOURCE__
38#define __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_SOURCE__
39
40#include "common/common.h"
54
55#ifdef UG_PARALLEL
58#endif
59
60namespace ug{
61namespace Electromagnetism{
62
64
121template <typename TDomain, typename TAlgebra>
123{
126
127public:
129 typedef TDomain domain_type;
130
132 typedef typename TDomain::grid_type grid_type;
133
135 typedef typename domain_type::subset_handler_type subset_handler_type;
136
138 typedef TAlgebra algebra_type;
140 typedef typename TAlgebra::vector_type vector_type;
142 typedef typename TAlgebra::matrix_type matrix_type;
143
152
154 static const int WDim = TDomain::dim;
155
157 typedef typename TDomain::position_type position_type;
158
159private:
160
166
167public:
168
171 (
172 const char * ssNames,
173 const char * posSsNames,
174 const char * cutSsNames,
175 SmartPtr<ApproximationSpace<TDomain> > vertApproxSpace,
177 );
178
180 void set
181 (
182 const char * fctNames,
183 number I
184 );
185
187 void compute
188 (
190 );
191
193 std::string subsets () {return m_allSsNames;};
194
200
201private:
202
204 template <typename TElem>
206 (
207 const DoFDistribution & edgeDD,
208 aa_edge_flag_type & in_source
209 );
212 {
216
218 (
219 this_type * pThis,
220 const DoFDistribution & edgeDD,
221 aa_edge_flag_type & in_source
222 )
223 : m_pThis (pThis), m_edgeDD (edgeDD), m_in_source (in_source) {}
224
225 template <typename TElem> void operator() (TElem &)
226 {
227 m_pThis->template mark_source_edges<TElem> (m_edgeDD, m_in_source);
228 }
229 };
230
233 (
234 pot_gf_type & pot_u
235 );
236
239 (
240 const DoFDistribution & vertDD,
241 pot_vector_type & src_pot,
242 DoFDistribution & edgeDD,
243 size_t func,
244 number value,
245 vector_type & src_field
246 );
247
249 template <typename TElem>
251 {
252 public:
254 static const int numCorners = ref_elem_type::numCorners;
255
257 static void stiffness
258 (
259 GridObject * elem,
260 const position_type vCornerCoords [],
261 number loc_A [numCorners] [numCorners]
262 );
263
265 static bool bnd_corners
266 (
267 GridObject * elem,
268 const SubsetGroup & bndGrp,
269 bool bnd_flag [numCorners]
270 );
271 };
272
274
279 class AuxLaplaceLocAss : public IElemDisc<TDomain>
280 {
281 private:
283 static const size_t _C_ = 0;
284
287
288 public:
291 (
292 NedelecLoopCurrent & master
293 );
294
296 virtual void prepare_setting
297 (
298 const std::vector<LFEID> & vLfeID,
299 bool bNonRegular
300 );
301
302 private:
305 template <typename TElem>
306 void prepare_element_loop (ReferenceObjectID roid, int si);
307
308 template <typename TElem>
309 void ass_JA_elem (LocalMatrix & J, const LocalVector & u, GridObject * elem, const position_type vCornerCoords []);
310
311 template <typename TElem>
312 void ass_rhs_elem (LocalVector & d, GridObject * elem, const position_type vCornerCoords []);
313
314 template <typename TElem>
315 void prepare_element (const LocalVector & u, GridObject * elem, ReferenceObjectID roid, const position_type vCornerCoords [])
316 {};
317 template <typename TElem>
319 {};
320 template <typename TElem>
321 void ass_JM_elem (LocalMatrix & J, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
322 {
323 UG_THROW ("NedelecProject: Attempt to use nonstationary discretization for the potential.");
324 };
325 template <typename TElem>
326 void ass_dA_elem (LocalVector & d, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
327 {
328 UG_THROW ("NedelecProject: Attempt to assemble non-linear system for the potential.");
329 };
330 template <typename TElem>
331 void ass_dM_elem (LocalVector & d, const LocalVector & u, GridObject * elem, const position_type vCornerCoords [])
332 {
333 UG_THROW ("NedelecProject: Attempt to use nonstationary discretization for the potential.");
334 };
336
340
344 template< typename TElem > void operator () (TElem&)
345 {m_pThis->register_loc_discr_func<TElem> ();}
346 };
347
348 template <typename TElem>
351
352 private:
353
356
361
364 };
365
367
373 class OutOfSource : public IDomainConstraint<TDomain, TPotAlgebra>
374 {
375 private:
378
381
382 private:
383
385 template <typename TElem>
387 (
388 const TDomain * dom
389 );
392 {
394 const TDomain * m_pDom;
395
397 (
398 OutOfSource * pThis,
399 const TDomain * pDom
400 )
401 : m_pThis (pThis), m_pDom (pDom) {}
402
403 template <typename TElem> void operator() (TElem &)
404 {
405 m_pThis->template mark_source_vertices_elem_type<TElem> (m_pDom);
406 }
407 };
410 (
411 const TDomain * dom
412 );
413
415 void adjust_matrix
416 (
417 const DoFDistribution & vertDD,
418 pot_matrix_type & A
419 );
420
422 void adjust_vector
423 (
424 const DoFDistribution & vertDD,
425 pot_vector_type & u
426 );
427
428 public:
431 (
432 const DoFDistribution & vertDD,
433 pot_vector_type & u
434 );
435
436 public:
439 (
440 NedelecLoopCurrent & master
441 );
442
444 void init
445 (
446 domain_type * dom,
447 a_vert_flag_type & a_in_source
448 )
449 {
450 MultiGrid * mg = dom->grid().get ();
451 m_in_source.access (*mg, a_in_source);
453# ifdef UG_PARALLEL
454 AttachmentAllReduce<Vertex> (*mg, a_in_source, PCL_RO_LOR);
455# endif
456 }
457
460 (
461 pot_matrix_type & J,
462 const pot_vector_type & u,
464 int type,
465 number time = 0.0,
467 const number s_a0 = 1.0
468 )
469 {
470 adjust_matrix (* dd.get (), J);
471 }
472
475 (
476 pot_vector_type & d,
477 const pot_vector_type & u,
479 int type,
480 number time = 0.0,
482 const std::vector<number> * vScaleMass = NULL,
483 const std::vector<number> * vScaleStiff = NULL
484 )
485 {
486 adjust_vector (* dd.get (), d);
487 }
488
491 (
492 pot_vector_type & u,
494 int type,
495 number time = 0.0
496 )
497 {
498 adjust_vector (* dd.get (), u);
499 }
500
503 (
504 pot_matrix_type & A,
505 pot_vector_type & b,
507 int type,
508 number time = 0.0
509 )
510 {
511 adjust_matrix (* dd.get (), A);
512 adjust_vector (* dd.get (), b);
513 }
514
517 (
518 pot_vector_type & b,
519 const pot_vector_type & u,
521 int type,
522 number time = 0.0
523 )
524 {
525 adjust_vector (* dd.get (), b);
526 }
527
529 int type () const {return CT_DIRICHLET;}
530
531 private:
532
535
538 };
539
541
546 class ZeroAverage : public ITransferPostProcess<TDomain, TPotAlgebra>
547 {
548 public:
549
552
555 {
556 m_spOoS->set_zero_average (* sp_gf->dd (), * sp_gf);
557 }
558
559 private:
560
563 };
564
566 template <typename TElem>
567 void get_flux_of_pot
568 (
569 const TDomain & domain,
570 const pot_vector_type & pot,
571 const DoFDistribution & vertDD,
572 number & flux
573 );
576 {
582
584 (
585 this_type * pThis,
586 const TDomain & domain,
587 const pot_vector_type & pot,
588 const DoFDistribution & vertDD,
589 number & flux
590 )
591 : m_pThis (pThis), m_domain (domain), m_pot (pot), m_vertDD (vertDD), m_flux (flux)
592 {
593 m_flux = 0;
594 }
595
596 template <typename TElem> void operator() (TElem &)
597 {
598 m_pThis->template get_flux_of_pot<TElem> (m_domain, m_pot, m_vertDD, m_flux);
599 }
600 };
601
602private:
603
605 struct TSrcData
606 {
607 std::string fctNames;
609
611 TSrcData (const char * the_fctNames, number the_I)
612 : fctNames (the_fctNames), I (the_I) {};
613 };
614
616 std::vector<TSrcData> m_vSrcData;
617
619 std::string m_allSsNames;
621 std::string m_posSsNames;
623 std::string m_cutSsNames;
624
631
634
642
653};
654
655} // end namespace Electromagnetism
656} // end namespace ug
657
658#include "nedelec_source_impl.h"
659
660#endif // __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_PROJECT__
661
662/* End of File */
Class for local assembling of the auxiliary Laplace operator.
Definition nedelec_source.h:280
bool m_in_pos_subset
whether we assemble in the subset specifying the positive direction
Definition nedelec_source.h:363
NedelecLoopCurrent & m_master
the master object of the projection
Definition nedelec_source.h:355
static const size_t _C_
a name for the component in the grid functions
Definition nedelec_source.h:283
SubsetGroup m_posSsGrp
group of the positive direction subsets
Definition nedelec_source.h:358
void ass_dM_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_source.h:331
void prepare_element(const LocalVector &u, GridObject *elem, ReferenceObjectID roid, const position_type vCornerCoords[])
Definition nedelec_source.h:315
void ass_JM_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_source.h:321
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_source_impl.h:496
void register_loc_discr_func()
Definition nedelec_source_impl.h:464
void register_all_loc_discr_funcs()
Definition nedelec_source_impl.h:452
SubsetGroup m_cutSsGrp
group of the surfaces at the cut
Definition nedelec_source.h:360
static const size_t maxCorners
maximum number of corners of an element
Definition nedelec_source.h:286
void finish_element_loop()
Definition nedelec_source.h:318
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegular)
switches between non-regular and regular grids
Definition nedelec_source_impl.h:432
void prepare_element_loop(ReferenceObjectID roid, int si)
Definition nedelec_source_impl.h:484
void ass_rhs_elem(LocalVector &d, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_source_impl.h:519
void ass_dA_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
Definition nedelec_source.h:326
Helper class for assembling the local stiffness matrix of Laplacian in various routines.
Definition nedelec_source.h:251
static void stiffness(GridObject *elem, const position_type vCornerCoords[], number loc_A[numCorners][numCorners])
Computes the local discretization of the Laplace operator.
Definition nedelec_source_impl.h:355
reference_element_traits< TElem >::reference_element_type ref_elem_type
Definition nedelec_source.h:253
static const int numCorners
Definition nedelec_source.h:254
static bool bnd_corners(GridObject *elem, const SubsetGroup &bndGrp, bool bnd_flag[numCorners])
Checks whether corners are on boundary (returns true if yes)
Definition nedelec_source_impl.h:389
Constraint that sets the problem to 0-identity out of the source.
Definition nedelec_source.h:374
DoFDistribution::traits< Vertex >::const_iterator t_vert_iterator
Iterator over vertices.
Definition nedelec_source.h:380
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_source.h:491
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_source.h:503
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_source.h:475
aa_vert_flag_type m_in_source
attachment accessor for the current vertex flags
Definition nedelec_source.h:537
int type() const
returns the type of the constraints
Definition nedelec_source.h:529
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_source.h:517
void adjust_vector(const DoFDistribution &vertDD, pot_vector_type &u)
sets to 0 all the entries of a vector that do not belong to the closure of the source subdomain
Definition nedelec_source_impl.h:650
void adjust_matrix(const DoFDistribution &vertDD, pot_matrix_type &A)
sets to identity all the matrix rows that do not belong to the closure of the source subdomain
Definition nedelec_source_impl.h:621
NedelecLoopCurrent & m_master
the master object of the projection
Definition nedelec_source.h:534
void init(domain_type *dom, a_vert_flag_type &a_in_source)
initializer
Definition nedelec_source.h:445
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_source.h:460
static const size_t maxCorners
maximum number of corners of an element
Definition nedelec_source.h:377
void mark_source_vertices_elem_type(const TDomain *dom)
marks the vertices belonging to elements in the source (for one element type)
Definition nedelec_source_impl.h:564
void mark_source_vertices(const TDomain *dom)
marks the vertices belonging to elements in the source (for all element types)
Definition nedelec_source_impl.h:600
void set_zero_average(const DoFDistribution &vertDD, pot_vector_type &u)
sets the arithmetic average of the solution vector to zero
Definition nedelec_source_impl.h:678
Postprocessor to eliminate the kernel parts in the solution.
Definition nedelec_source.h:547
virtual void post_process(SmartPtr< pot_gf_type > sp_gf)
Sets the arithmetical average to zero.
Definition nedelec_source.h:554
SmartPtr< OutOfSource > m_spOoS
access to the inner vertex marks
Definition nedelec_source.h:562
ZeroAverage(SmartPtr< OutOfSource > spOoS)
Constructor.
Definition nedelec_source.h:551
Class for computation of loop currents.
Definition nedelec_source.h:123
Grid::VertexAttachmentAccessor< a_vert_flag_type > aa_vert_flag_type
Definition nedelec_source.h:165
ABool a_vert_flag_type
Definition nedelec_source.h:164
SmartPtr< AssembledLinearOperator< TPotAlgebra > > m_auxLaplaceOp
Matrix of the discretization.
Definition nedelec_source.h:650
void compute(SmartPtr< GridFunction< TDomain, TAlgebra > > sp_u)
Computation of the source (only if 'init'ialized): distributes the potential.
Definition nedelec_source_impl.h:121
std::string subsets()
Returns the subsets where the source is defined.
Definition nedelec_source.h:193
SmartPtr< OutOfSource > m_outOfSource
Extension of the matrices and vectors to the whole domain.
Definition nedelec_source.h:638
SubsetGroup m_posSsGrp
Group of the subsets of the positive direction.
Definition nedelec_source.h:628
SubsetGroup m_cutSsGrp
Group of the surfaces of the cut of the loop.
Definition nedelec_source.h:630
void distribute_source_potential(const DoFDistribution &vertDD, pot_vector_type &src_pot, DoFDistribution &edgeDD, size_t func, number value, vector_type &src_field)
Computes the gradient of the potential.
Definition nedelec_source_impl.h:281
std::string m_allSsNames
Names of all the subsets of the source.
Definition nedelec_source.h:619
std::string m_cutSsNames
Names of the surfaces of the cut of the loop.
Definition nedelec_source.h:623
TDomain::grid_type grid_type
Type of Grid:
Definition nedelec_source.h:132
TAlgebra::vector_type vector_type
Type of Vector (for the Nedelec-element-based grid functions)
Definition nedelec_source.h:140
SmartPtr< ITransferPostProcess< TDomain, TPotAlgebra > > zero_average()
Returns the postprocess routine to eliminate the kernel part of the potential.
Definition nedelec_source.h:196
AChar a_edge_flag_type
Type of the attachment and its accessor for flags.
Definition nedelec_source.h:162
void compute_potential(pot_gf_type &pot_u)
Computes the potential of the source.
Definition nedelec_source_impl.h:246
std::vector< TSrcData > m_vSrcData
Array of the electric current data.
Definition nedelec_source.h:616
TPotAlgebra::matrix_type pot_matrix_type
Matrix type for the potential space.
Definition nedelec_source.h:149
domain_type::subset_handler_type subset_handler_type
Type of subset handler.
Definition nedelec_source.h:135
TDomain domain_type
Type of Domain.
Definition nedelec_source.h:129
TDomain::position_type position_type
position type
Definition nedelec_source.h:157
CPUAlgebra TPotAlgebra
The auxiliary algebra type for the space of the potential. (Note: It should be scalar....
Definition nedelec_source.h:145
std::string m_posSsNames
Names of the subsets of the positive direction.
Definition nedelec_source.h:621
SmartPtr< ZeroAverage > m_zeroAverage
Elimination of the kernal parts in the solution.
Definition nedelec_source.h:640
TAlgebra algebra_type
Type of algebra (for the Nedelec-element-based grid functions)
Definition nedelec_source.h:138
SmartPtr< AuxLaplaceLocAss > m_auxLocLaplace
Local discretization of the auxiliary equations.
Definition nedelec_source.h:636
TAlgebra::matrix_type matrix_type
Type of Vector (for the Nedelec-element-based grid functions)
Definition nedelec_source.h:142
GridFunction< TDomain, TPotAlgebra > pot_gf_type
Grid function type for the potential space.
Definition nedelec_source.h:151
void set(const char *fctNames, number I)
Setting the electric current value.
Definition nedelec_source_impl.h:108
SmartPtr< DomainDiscretization< TDomain, TPotAlgebra > > m_auxLaplaceAss
Global discretization of the auxiliary equations.
Definition nedelec_source.h:648
NedelecLoopCurrent< TDomain, TAlgebra > this_type
This type.
Definition nedelec_source.h:125
void get_flux_of_pot(const TDomain &domain, const pot_vector_type &pot, const DoFDistribution &vertDD, number &flux)
Computes the flux of of the gradient of the potential over the cut (for one type of elements)
Definition nedelec_source_impl.h:733
void mark_source_edges(const DoFDistribution &edgeDD, aa_edge_flag_type &in_source)
Marks edges that belong to the loop source domain (for one type of elements)
Definition nedelec_source_impl.h:199
static const int WDim
world dimention
Definition nedelec_source.h:154
Grid::EdgeAttachmentAccessor< a_edge_flag_type > aa_edge_flag_type
Definition nedelec_source.h:163
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
Vertex-centered approximation space.
Definition nedelec_source.h:633
SubsetGroup m_allSsGrp
Group of all the subsets of the source.
Definition nedelec_source.h:626
TPotAlgebra::vector_type pot_vector_type
Vector type for the potential space.
Definition nedelec_source.h:147
SmartPtr< ILinearOperatorInverse< pot_vector_type > > m_potSolver
Solver for the auxliliary equations.
Definition nedelec_source.h:652
bool access(Grid &grid, TAttachment &a)
ConstSmartPtr< DoFDistribution > dd(const GridLevel &gl) const
#define PCL_RO_LOR
#define UG_THROW(msg)
double number
CT_DIRICHLET
ReferenceObjectID
RegisterLocalDiscr(AuxLaplaceLocAss *pThis)
Definition nedelec_source.h:342
Helper class for computation of the flux of potential over the cut.
Definition nedelec_source.h:576
const domain_type & m_domain
Definition nedelec_source.h:578
GetFluxOfPotential(this_type *pThis, const TDomain &domain, const pot_vector_type &pot, const DoFDistribution &vertDD, number &flux)
Definition nedelec_source.h:584
number & m_flux
Definition nedelec_source.h:581
const DoFDistribution & m_vertDD
Definition nedelec_source.h:580
const pot_vector_type & m_pot
Definition nedelec_source.h:579
void operator()(TElem &)
Definition nedelec_source.h:596
this_type * m_pThis
Definition nedelec_source.h:577
Helper class for marking the edges in the source.
Definition nedelec_source.h:212
MarkSourceEdges(this_type *pThis, const DoFDistribution &edgeDD, aa_edge_flag_type &in_source)
Definition nedelec_source.h:218
void operator()(TElem &)
Definition nedelec_source.h:225
this_type * m_pThis
Definition nedelec_source.h:213
const DoFDistribution & m_edgeDD
Definition nedelec_source.h:214
aa_edge_flag_type & m_in_source
Definition nedelec_source.h:215
Helper class for marking the vertices in the closure of the source subdomain.
Definition nedelec_source.h:392
MarkSourceVertices(OutOfSource *pThis, const TDomain *pDom)
Definition nedelec_source.h:397
Structure for keeping electric current data.
Definition nedelec_source.h:606
TSrcData(const char *the_fctNames, number the_I)
Constructor.
Definition nedelec_source.h:611
number I
value of the electric current
Definition nedelec_source.h:608
std::string fctNames
function names for the value
Definition nedelec_source.h:607