dune-pdelab  2.5-dev
maxwelldg.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil -*-
2 #ifndef DUNE_PDELAB_LOCALOPERATOR_MAXWELLDG_HH
3 #define DUNE_PDELAB_LOCALOPERATOR_MAXWELLDG_HH
4 
5 #include<vector>
6 
7 #include<dune/common/exceptions.hh>
8 #include<dune/common/fvector.hh>
9 
10 #include <dune/geometry/referenceelements.hh>
11 
21 
22 #include"maxwellparameter.hh"
23 
24 namespace Dune {
25  namespace PDELab {
26 
27 
28  template<int dim>
30  {};
31 
37  template<>
39  {
40  enum { dim = 3 };
41  public:
42 
52  template<typename T1, typename T2, typename T3>
53  static void eigenvalues (T1 eps, T1 mu, const Dune::FieldVector<T2,2*dim>& e)
54  {
55  T1 s = 1.0/sqrt(mu*eps); //speed of light s = 1/sqrt(\mu \epsilon)
56  e[0] = s;
57  e[1] = s;
58  e[2] = -s;
59  e[3] = -s;
60  e[4] = 0;
61  e[5] = 0;
62  }
63 
75  template<typename T1, typename T2, typename T3>
76  static void eigenvectors (T1 eps, T1 mu, const Dune::FieldVector<T2,dim>& n, Dune::FieldMatrix<T3,2*dim,2*dim>& R)
77  {
78  T1 a=n[0], b=n[1], c=n[2];
79 
80  Dune::FieldVector<T2,dim> re, im;
81  if (std::abs(c)<0.5)
82  {
83  re[0]=a*c; re[1]=b*c; re[2]=c*c-1;
84  im[0]=-b; im[1]=a; im[2]=0;
85  }
86  else
87  {
88  re[0]=a*b; re[1]=b*b-1; re[2]=b*c;
89  im[0]=c; im[1]=0.0; im[2]=-a;
90  }
91 
92  // \lambda_0,1 = s
93  R[0][0] = re[0]; R[0][1] = -im[0];
94  R[1][0] = re[1]; R[1][1] = -im[1];
95  R[2][0] = re[2]; R[2][1] = -im[2];
96  R[3][0] = im[0]; R[3][1] = re[0];
97  R[4][0] = im[1]; R[4][1] = re[1];
98  R[5][0] = im[2]; R[5][1] = re[2];
99 
100  // \lambda_2,3 = -s
101  R[0][2] = im[0]; R[0][3] = re[0];
102  R[1][2] = im[1]; R[1][3] = re[1];
103  R[2][2] = im[2]; R[2][3] = re[2];
104  R[3][2] = re[0]; R[3][3] = -im[0];
105  R[4][2] = re[1]; R[4][3] = -im[1];
106  R[5][2] = re[2]; R[5][3] = -im[2];
107 
108  // \lambda_4,5 = 0
109  R[0][4] = a; R[0][5] = 0;
110  R[1][4] = b; R[1][5] = 0;
111  R[2][4] = c; R[2][5] = 0;
112  R[3][4] = 0; R[3][5] = a;
113  R[4][4] = 0; R[4][5] = b;
114  R[5][4] = 0; R[5][5] = c;
115 
116  // apply scaling
117  T1 weps=sqrt(eps);
118  T1 wmu=sqrt(mu);
119  for (std::size_t i=0; i<3; i++)
120  for (std::size_t j=0; j<6; j++)
121  R[i][j] *= weps;
122  for (std::size_t i=3; i<6; i++)
123  for (std::size_t j=0; j<6; j++)
124  R[i][j] *= wmu;
125 
126  return;
127 
128  // if (std::abs(std::abs(c)-1)<1e-10)
129  // {
130  // if (c>0)
131  // {
132  // // \lambda_0,1 = s
133  // R[0][0] = 0; R[0][1] = 1;
134  // R[1][0] = -1; R[1][1] = 0;
135  // R[2][0] = 0; R[2][1] = 0;
136  // R[3][0] = 1; R[3][1] = 0;
137  // R[4][0] = 0; R[4][1] = 1;
138  // R[5][0] = 0; R[5][1] = 0;
139 
140  // // \lambda_2,3 = -s
141  // R[0][2] = -1; R[0][3] = 0;
142  // R[1][2] = 0; R[1][3] = 1;
143  // R[2][2] = 0; R[2][3] = 0;
144  // R[3][2] = 0; R[3][3] = 1;
145  // R[4][2] = 1; R[4][3] = 0;
146  // R[5][2] = 0; R[5][3] = 0;
147  // }
148  // else
149  // {
150  // // \lambda_0,1 = s
151  // R[0][0] = -1; R[0][1] = 0;
152  // R[1][0] = 0; R[1][1] = 1;
153  // R[2][0] = 0; R[2][1] = 0;
154  // R[3][0] = 0; R[3][1] = 1;
155  // R[4][0] = 1; R[4][1] = 0;
156  // R[5][0] = 0; R[5][1] = 0;
157 
158  // // \lambda_2,3 = -s
159  // R[0][2] = 0; R[0][3] = 1;
160  // R[1][2] = -1; R[1][3] = 0;
161  // R[2][2] = 0; R[2][3] = 0;
162  // R[3][2] = 1; R[3][3] = 0;
163  // R[4][2] = 0; R[4][3] = 1;
164  // R[5][2] = 0; R[5][3] = 0;
165  // }
166  // }
167  // else if (std::abs(std::abs(b)-1)<1e-10)
168  // {
169  // if (b>0)
170  // {
171  // // \lambda_0,1 = s
172  // R[0][0] = -1; R[0][1] = 0;
173  // R[1][0] = 0; R[1][1] = 0;
174  // R[2][0] = 0; R[2][1] = 1;
175  // R[3][0] = 0; R[3][1] = 1;
176  // R[4][0] = 0; R[4][1] = 0;
177  // R[5][0] = 1; R[5][1] = 0;
178 
179  // // \lambda_2,3 = -s
180  // R[0][2] = 0; R[0][3] = 1;
181  // R[1][2] = 0; R[1][3] = 0;
182  // R[2][2] = -1; R[2][3] = 0;
183  // R[3][2] = 1; R[3][3] = 0;
184  // R[4][2] = 0; R[4][3] = 0;
185  // R[5][2] = 0; R[5][3] = 1;
186  // }
187  // else
188  // {
189  // // \lambda_0,1 = s
190  // R[0][0] = 0; R[0][1] = 1;
191  // R[1][0] = 0; R[1][1] = 0;
192  // R[2][0] = -1; R[2][1] = 0;
193  // R[3][0] = 1; R[3][1] = 0;
194  // R[4][0] = 0; R[4][1] = 0;
195  // R[5][0] = 0; R[5][1] = 1;
196 
197  // // \lambda_2,3 = -s
198  // R[0][2] = -1; R[0][3] = 0;
199  // R[1][2] = 0; R[1][3] = 0;
200  // R[2][2] = 0; R[2][3] = 1;
201  // R[3][2] = 0; R[3][3] = 1;
202  // R[4][2] = 0; R[4][3] = 0;
203  // R[5][2] = 1; R[5][3] = 0;
204  // }
205  // }
206  // else if (std::abs(std::abs(a)-1)<1e-10)
207  // {
208  // if (a>0)
209  // {
210  // // \lambda_0,1 = s
211  // R[0][0] = 0; R[0][1] = 0;
212  // R[1][0] = 0; R[1][1] = 1;
213  // R[2][0] = -1; R[2][1] = 0;
214  // R[3][0] = 0; R[3][1] = 0;
215  // R[4][0] = 1; R[4][1] = 0;
216  // R[5][0] = 0; R[5][1] = 1;
217 
218  // // \lambda_2,3 = -s
219  // R[0][2] = 0; R[0][3] = 0;
220  // R[1][2] = -1; R[1][3] = 0;
221  // R[2][2] = 0; R[2][3] = 1;
222  // R[3][2] = 0; R[3][3] = 0;
223  // R[4][2] = 0; R[4][3] = 1;
224  // R[5][2] = 1; R[5][3] = 0;
225  // }
226  // else
227  // {
228  // // \lambda_0,1 = s
229  // R[0][0] = 0; R[0][1] = 0;
230  // R[1][0] = -1; R[1][1] = 0;
231  // R[2][0] = 0; R[2][1] = 1;
232  // R[3][0] = 0; R[3][1] = 0;
233  // R[4][0] = 0; R[4][1] = 1;
234  // R[5][0] = 1; R[5][1] = 0;
235 
236  // // \lambda_2,3 = -s
237  // R[0][2] = 0; R[0][3] = 0;
238  // R[1][2] = 0; R[1][3] = 1;
239  // R[2][2] = -1; R[2][3] = 0;
240  // R[3][2] = 0; R[3][3] = 0;
241  // R[4][2] = 1; R[4][3] = 0;
242  // R[5][2] = 0; R[5][3] = 1;
243  // }
244  // }
245  // else
246  // {
247  // DUNE_THROW(Dune::Exception,"need axiparallel grids for now!");
248 
249  // // \lambda_0,1 = s
250  // R[0][0] = b; R[0][1] = -(b*b+c*c);
251  // R[1][0] = -a; R[1][1] = a*b;
252  // R[2][0] = 0; R[2][1] = a*c;
253  // R[3][0] = a*c; R[3][1] = 0;
254  // R[4][0] = b*c; R[4][1] = -c;
255  // R[5][0] = -(a*a+b*b); R[5][1] = b;
256 
257  // // \lambda_2,3 = -s
258  // R[0][2] = -b; R[0][3] = -(b*b+c*c);
259  // R[1][2] = a; R[1][3] = a*b;
260  // R[2][2] = 0; R[2][3] = a*c;
261  // R[3][2] = a*c; R[3][3] = 0;
262  // R[4][2] = b*c; R[4][3] = c;
263  // R[5][2] = -(a*a+b*b); R[5][3] = -b;
264  // }
265 
266  // // \lambda_4,5 = 0
267  // R[0][4] = 0; R[0][5] = a;
268  // R[1][4] = 0; R[1][5] = b;
269  // R[2][4] = 0; R[2][5] = c;
270  // R[3][4] = a; R[3][5] = 0;
271  // R[4][4] = b; R[4][5] = 0;
272  // R[5][4] = c; R[5][5] = 0;
273 
274  // // apply scaling
275  // T1 weps=sqrt(eps);
276  // T1 wmu=sqrt(mu);
277  // for (std::size_t i=0; i<3; i++)
278  // for (std::size_t j=0; j<6; j++)
279  // R[i][j] *= weps;
280  // for (std::size_t i=3; i<6; i++)
281  // for (std::size_t j=0; j<6; j++)
282  // R[i][j] *= wmu;
283 
284  //std::cout << R << std::endl;
285 
286  }
287  };
288 
303  template<typename T, typename FEM>
305  public NumericalJacobianApplyVolume<DGMaxwellSpatialOperator<T,FEM> >,
306  public NumericalJacobianVolume<DGMaxwellSpatialOperator<T,FEM> >,
307  public NumericalJacobianApplySkeleton<DGMaxwellSpatialOperator<T,FEM> >,
308  public NumericalJacobianSkeleton<DGMaxwellSpatialOperator<T,FEM> >,
309  public NumericalJacobianApplyBoundary<DGMaxwellSpatialOperator<T,FEM> >,
310  public NumericalJacobianBoundary<DGMaxwellSpatialOperator<T,FEM> >,
311  public FullSkeletonPattern,
312  public FullVolumePattern,
314  public InstationaryLocalOperatorDefaultMethods<typename T::Traits::RangeFieldType>
315  {
316  enum { dim = T::Traits::GridViewType::dimension };
317 
318  public:
319  // pattern assembly flags
320  enum { doPatternVolume = true };
321  enum { doPatternSkeleton = true };
322 
323  // residual assembly flags
324  enum { doAlphaVolume = true };
325  enum { doAlphaSkeleton = true };
326  enum { doAlphaBoundary = true };
327  enum { doLambdaVolume = true };
328 
329  // ! constructor
330  DGMaxwellSpatialOperator (T& param_, int overintegration_=0)
331  : param(param_), overintegration(overintegration_), cache(20)
332  {
333  }
334 
335  // volume integral depending on test and ansatz functions
336  template<typename EG, typename LFSU, typename X, typename LFSV, typename R>
337  void alpha_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv, R& r) const
338  {
339  // Define types
340  using namespace Indices;
341  using DGSpace = TypeTree::Child<LFSV,_0>;
342  using RF = typename DGSpace::Traits::FiniteElementType::Traits
343  ::LocalBasisType::Traits::RangeFieldType;
344  using size_type = typename DGSpace::Traits::SizeType;
345 
346  // paranoia check number of number of components
347  static_assert(TypeTree::StaticDegree<LFSV>::value == dim*2,
348  "need exactly dim*2 components!");
349 
350  // get local function space that is identical for all components
351  const auto& dgspace = child(lfsv,_0);
352 
353  // Reference to cell
354  const auto& cell = eg.entity();
355 
356  // Get geometry
357  auto geo = eg.geometry();
358 
359  // evaluate parameters (assumed constant per element)
360  auto ref_el = referenceElement(geo);
361  auto localcenter = ref_el.position(0,0);
362  auto mu = param.mu(cell,localcenter);
363  auto eps = param.eps(cell,localcenter);
364  auto sigma = param.sigma(cell,localcenter);
365  auto muinv = 1.0/mu;
366  auto epsinv = 1.0/eps;
367 
368  //std::cout << "alpha_volume center=" << eg.geometry().center() << std::endl;
369 
370  // Transformation
371  typename EG::Geometry::JacobianInverseTransposed jac;
372 
373  // Initialize vectors outside for loop
374  Dune::FieldVector<RF,dim*2> u(0.0);
375  std::vector<Dune::FieldVector<RF,dim> > gradphi(dgspace.size());
376 
377  // loop over quadrature points
378  const int order = dgspace.finiteElement().localBasis().order();
379  const int intorder = overintegration+2*order;
380  for (const auto &qp : quadratureRule(geo,intorder))
381  {
382  // evaluate basis functions
383  const auto& phi = cache[order].evaluateFunction
384  (qp.position(), dgspace.finiteElement().localBasis());
385 
386  // evaluate state vector u
387  for (size_type k=0; k<dim*2; k++){ // for all components
388  u[k] = 0.0;
389  for (size_type j=0; j<dgspace.size(); j++) // for all basis functions
390  u[k] += x(lfsv.child(k),j)*phi[j];
391  }
392  //std::cout << " u at " << qp.position() << " : " << u << std::endl;
393 
394  // evaluate gradient of basis functions (we assume Galerkin method lfsu=lfsv)
395  const auto& js = cache[order].evaluateJacobian
396  (qp.position(), dgspace.finiteElement().localBasis());
397 
398  // compute global gradients
399  jac = geo.jacobianInverseTransposed(qp.position());
400  for (size_type i=0; i<dgspace.size(); i++)
401  jac.mv(js[i][0],gradphi[i]);
402 
403  // integrate
404  auto factor = qp.weight() * geo.integrationElement(qp.position());
405 
406  Dune::FieldMatrix<RF,dim*2,dim> F;
407  F[0][0] = 0; F[0][1] = -muinv*u[5]; F[0][2] = muinv*u[4];
408  F[1][0] = muinv*u[5]; F[1][1] = 0; F[1][2] = -muinv*u[3];
409  F[2][0] =-muinv*u[4]; F[2][1] = muinv*u[3]; F[2][2] = 0;
410  F[3][0] = 0; F[3][1] = epsinv*u[2]; F[3][2] = -epsinv*u[1];
411  F[4][0] = -epsinv*u[2]; F[4][1] = 0; F[4][2] = epsinv*u[0];
412  F[5][0] = epsinv*u[1]; F[5][1] = -epsinv*u[0]; F[5][2] = 0;
413 
414  // for all components of the system
415  for (size_type i=0; i<dim*2; i++)
416  // for all test functions of this component
417  for (size_type k=0; k<dgspace.size(); k++)
418  // for all dimensions
419  for (size_type j=0; j<dim; j++)
420  r.accumulate(lfsv.child(i),k,-F[i][j]*gradphi[k][j]*factor);
421 
422  // for the first half of the system
423  for (size_type i=0; i<dim; i++)
424  // for all test functions of this component
425  for (size_type k=0; k<dgspace.size(); k++)
426  r.accumulate(lfsv.child(i),k,(sigma/eps)*u[i]*phi[k]*factor);
427 
428  // std::cout << " residual: ";
429  // for (size_type i=0; i<r.size(); i++) std::cout << r[i] << " ";
430  // std::cout << std::endl;
431  }
432  }
433 
434  // skeleton integral depending on test and ansatz functions
435  // each face is only visited ONCE!
436  template<typename IG, typename LFSU, typename X, typename LFSV, typename R>
437  void alpha_skeleton (const IG& ig,
438  const LFSU& lfsu_s, const X& x_s, const LFSV& lfsv_s,
439  const LFSU& lfsu_n, const X& x_n, const LFSV& lfsv_n,
440  R& r_s, R& r_n) const
441  {
442  using std::max;
443  using std::sqrt;
444 
445  // Define types
446  using namespace Indices;
447  using DGSpace = TypeTree::Child<LFSV,_0>;
448  using DF = typename DGSpace::Traits::FiniteElementType::
449  Traits::LocalBasisType::Traits::DomainFieldType;
450  using RF = typename DGSpace::Traits::FiniteElementType::
451  Traits::LocalBasisType::Traits::RangeFieldType;
452  using size_type = typename DGSpace::Traits::SizeType;
453 
454  // get local function space that is identical for all components
455  const auto& dgspace_s = child(lfsv_s,_0);
456  const auto& dgspace_n = child(lfsv_n,_0);
457 
458  // References to inside and outside cells
459  const auto& cell_inside = ig.inside();
460  const auto& cell_outside = ig.outside();
461 
462  // Get geometries
463  auto geo = ig.geometry();
464  auto geo_inside = cell_inside.geometry();
465  auto geo_outside = cell_outside.geometry();
466 
467  // Geometry of intersection in local coordinates of inside_cell and outside_cell
468  auto geo_in_inside = ig.geometryInInside();
469  auto geo_in_outside = ig.geometryInOutside();
470 
471  // evaluate speed of sound (assumed constant per element)
472  auto ref_el_inside = referenceElement(geo_inside);
473  auto ref_el_outside = referenceElement(geo_outside);
474  auto local_inside = ref_el_inside.position(0,0);
475  auto local_outside = ref_el_outside.position(0,0);
476  // This is wrong -- this approach (with A- and A+) does not allow
477  // position-dependent eps and mu, so we should not allow the parameter
478  // class to specify them in a position-dependent manner. See my
479  // (Jorrit Fahlke) dissertation on how to do it right.
480  auto mu_s = param.mu(cell_inside,local_inside);
481  auto mu_n = param.mu(cell_outside,local_outside);
482  auto eps_s = param.eps(cell_inside,local_inside);
483  auto eps_n = param.eps(cell_outside,local_outside);
484  //auto sigma_s = param.sigma(ig.inside(),local_inside);
485  //auto sigma_n = param.sigma(ig.outside(),local_outside);
486 
487  // normal: assume faces are planar
488  const auto& n_F = ig.centerUnitOuterNormal();
489 
490  // compute A+ (outgoing waves)
491  Dune::FieldMatrix<DF,dim*2,dim*2> R_s;
492  MaxwellEigenvectors<dim>::eigenvectors(eps_s,mu_s,n_F,R_s);
493  Dune::FieldMatrix<DF,dim*2,dim*2> Dplus_s(0.0);
494  Dplus_s[0][0] = 1.0/sqrt(eps_s*mu_s);
495  Dplus_s[1][1] = 1.0/sqrt(eps_s*mu_s);
496  Dune::FieldMatrix<DF,dim*2,dim*2> Aplus_s(R_s);
497  Aplus_s.rightmultiply(Dplus_s);
498  R_s.invert();
499  Aplus_s.rightmultiply(R_s);
500 
501  // compute A- (incoming waves)
502  Dune::FieldMatrix<DF,dim*2,dim*2> R_n;
503  MaxwellEigenvectors<dim>::eigenvectors(eps_n,mu_n,n_F,R_n);
504  Dune::FieldMatrix<DF,dim*2,dim*2> Dminus_n(0.0);
505  Dminus_n[2][2] = -1.0/sqrt(eps_n*mu_n);
506  Dminus_n[3][3] = -1.0/sqrt(eps_n*mu_n);
507  Dune::FieldMatrix<DF,dim*2,dim*2> Aminus_n(R_n);
508  Aminus_n.rightmultiply(Dminus_n);
509  R_n.invert();
510  Aminus_n.rightmultiply(R_n);
511 
512  // Initialize vectors outside for loop
513  Dune::FieldVector<RF,dim*2> u_s(0.0);
514  Dune::FieldVector<RF,dim*2> u_n(0.0);
515  Dune::FieldVector<RF,dim*2> f(0.0);
516 
517  // std::cout << "alpha_skeleton center=" << ig.geometry().center() << std::endl;
518 
519  // loop over quadrature points
520  const int order_s = dgspace_s.finiteElement().localBasis().order();
521  const int order_n = dgspace_n.finiteElement().localBasis().order();
522  const int intorder = overintegration+1+2*max(order_s,order_n);
523  for (const auto& qp : quadratureRule(geo,intorder))
524  {
525  // position of quadrature point in local coordinates of elements
526  const auto& iplocal_s = geo_in_inside.global(qp.position());
527  const auto& iplocal_n = geo_in_outside.global(qp.position());
528 
529  // evaluate basis functions
530  const auto& phi_s = cache[order_s].evaluateFunction(iplocal_s,dgspace_s.finiteElement().localBasis());
531  const auto& phi_n = cache[order_n].evaluateFunction(iplocal_n,dgspace_n.finiteElement().localBasis());
532 
533  // evaluate u from inside and outside
534  for (size_type i=0; i<dim*2; i++){ // for all components
535  u_s[i] = 0.0;
536  for (size_type k=0; k<dgspace_s.size(); k++) // for all basis functions
537  u_s[i] += x_s(lfsv_s.child(i),k)*phi_s[k];
538  }
539  for (size_type i=0; i<dim*2; i++){ // for all components
540  u_n[i] = 0.0;
541  for (size_type k=0; k<dgspace_n.size(); k++) // for all basis functions
542  u_n[i] += x_n(lfsv_n.child(i),k)*phi_n[k];
543  }
544 
545  // compute numerical flux at integration point
546  f = 0.0;
547  Aplus_s.umv(u_s,f);
548  // std::cout << " after A_plus*u_s " << f << std::endl;
549  Aminus_n.umv(u_n,f);
550  // std::cout << " after A_minus*u_n " << f << std::endl;
551 
552  //std::cout << "f=" << f << " u_s=" << u_s << " u_n=" << u_n << std::endl;
553 
554  // integrate
555  auto factor = qp.weight() * geo.integrationElement(qp.position());
556  for (size_type k=0; k<dgspace_s.size(); k++) // loop over all vector-valued (!) basis functions (with identical components)
557  for (size_type i=0; i<dim*2; i++) // loop over all components
558  r_s.accumulate(lfsv_s.child(i),k,f[i]*phi_s[k]*factor);
559  for (size_type k=0; k<dgspace_n.size(); k++) // loop over all vector-valued (!) basis functions (with identical components)
560  for (size_type i=0; i<dim*2; i++) // loop over all components
561  r_n.accumulate(lfsv_n.child(i),k,-f[i]*phi_n[k]*factor);
562  }
563 
564  // std::cout << " residual_s: ";
565  // for (auto v : r_s) std::cout << v << " ";
566  // std::cout << std::endl;
567  // std::cout << " residual_n: ";
568  // for (auto v : r_n) std::cout << v << " ";
569  // std::cout << std::endl;
570  }
571 
572  // skeleton integral depending on test and ansatz functions
573  template<typename IG, typename LFSU, typename X, typename LFSV, typename R>
574  void alpha_boundary (const IG& ig,
575  const LFSU& lfsu_s, const X& x_s, const LFSV& lfsv_s,
576  R& r_s) const
577  {
578  // Define types
579  using namespace Indices;
580  using DGSpace = TypeTree::Child<LFSV,_0>;
581  using DF = typename DGSpace::Traits::FiniteElementType::
582  Traits::LocalBasisType::Traits::DomainFieldType;
583  using RF = typename DGSpace::Traits::FiniteElementType::
584  Traits::LocalBasisType::Traits::RangeFieldType;
585  using size_type = typename DGSpace::Traits::SizeType;
586 
587  // get local function space that is identical for all components
588  const auto& dgspace_s = child(lfsv_s,_0);
589 
590  // References to inside cell
591  const auto& cell_inside = ig.inside();
592 
593  // Get geometries
594  auto geo = ig.geometry();
595  auto geo_inside = cell_inside.geometry();
596 
597  // Geometry of intersection in local coordinates of inside_cell
598  auto geo_in_inside = ig.geometryInInside();
599 
600  // evaluate speed of sound (assumed constant per element)
601  auto ref_el_inside = referenceElement(geo_inside);
602  auto local_inside = ref_el_inside.position(0,0);
603  auto mu_s = param.mu(cell_inside,local_inside);
604  auto eps_s = param.eps(cell_inside,local_inside);
605  //auto sigma_s = param.sigma(ig.inside(),local_inside );
606 
607  // normal: assume faces are planar
608  const auto& n_F = ig.centerUnitOuterNormal();
609 
610  // compute A+ (outgoing waves)
611  Dune::FieldMatrix<DF,dim*2,dim*2> R_s;
612  MaxwellEigenvectors<dim>::eigenvectors(eps_s,mu_s,n_F,R_s);
613  Dune::FieldMatrix<DF,dim*2,dim*2> Dplus_s(0.0);
614  Dplus_s[0][0] = 1.0/sqrt(eps_s*mu_s);
615  Dplus_s[1][1] = 1.0/sqrt(eps_s*mu_s);
616  Dune::FieldMatrix<DF,dim*2,dim*2> Aplus_s(R_s);
617  Aplus_s.rightmultiply(Dplus_s);
618  R_s.invert();
619  Aplus_s.rightmultiply(R_s);
620 
621  // compute A- (incoming waves)
622  Dune::FieldMatrix<DF,dim*2,dim*2> R_n;
623  MaxwellEigenvectors<dim>::eigenvectors(eps_s,mu_s,n_F,R_n);
624  Dune::FieldMatrix<DF,dim*2,dim*2> Dminus_n(0.0);
625  Dminus_n[2][2] = -1.0/sqrt(eps_s*mu_s);
626  Dminus_n[3][3] = -1.0/sqrt(eps_s*mu_s);
627  Dune::FieldMatrix<DF,dim*2,dim*2> Aminus_n(R_n);
628  Aminus_n.rightmultiply(Dminus_n);
629  R_n.invert();
630  Aminus_n.rightmultiply(R_n);
631 
632  // Initialize vectors outside for loop
633  Dune::FieldVector<RF,dim*2> u_s(0.0);
634  Dune::FieldVector<RF,dim*2> u_n;
635  Dune::FieldVector<RF,dim*2> f(0.0);
636 
637  // std::cout << "alpha_boundary center=" << ig.geometry().center() << std::endl;
638 
639  // loop over quadrature points
640  const int order_s = dgspace_s.finiteElement().localBasis().order();
641  const int intorder = overintegration+1+2*order_s;
642  for(const auto &qp : quadratureRule(geo,intorder))
643  {
644  // position of quadrature point in local coordinates of elements
645  const auto& iplocal_s = geo_in_inside.global(qp.position());
646 
647  // evaluate basis functions
648  const auto& phi_s = cache[order_s].evaluateFunction
649  (iplocal_s,dgspace_s.finiteElement().localBasis());
650 
651  // evaluate u from inside and outside
652  for (size_type i=0; i<dim*2; i++){ // for all components
653  u_s[i] = 0.0;
654  for (size_type k=0; k<dgspace_s.size(); k++) // for all basis functions
655  u_s[i] += x_s(lfsv_s.child(i),k)*phi_s[k];
656  }
657  // std::cout << " u_s " << u_s << std::endl;
658 
659  // evaluate boundary condition
660  u_n = (param.g(ig.intersection(),qp.position(),u_s));
661  // std::cout << " u_n " << u_n << " bc: " << param.g(ig.intersection(),qp.position(),u_s) << std::endl;
662 
663  // compute numerical flux at integration point
664  f = 0.0;
665  Aplus_s.umv(u_s,f);
666  // std::cout << " after A_plus*u_s " << f << std::endl;
667  Aminus_n.umv(u_n,f);
668  // std::cout << " after A_minus*u_n " << f << std::endl;
669 
670  // integrate
671  auto factor = qp.weight() * geo.integrationElement(qp.position());
672  for (size_type k=0; k<dgspace_s.size(); k++) // loop over all vector-valued (!) basis functions (with identical components)
673  for (size_type i=0; i<dim*2; i++) // loop over all components
674  r_s.accumulate(lfsv_s.child(i),k,f[i]*phi_s[k]*factor);
675  }
676 
677  // std::cout << " residual_s: ";
678  // for (size_type i=0; i<r_s.size(); i++) std::cout << r_s[i] << " ";
679  // std::cout << std::endl;
680  }
681 
682  // volume integral depending only on test functions
683  template<typename EG, typename LFSV, typename R>
684  void lambda_volume (const EG& eg, const LFSV& lfsv, R& r) const
685  {
686  // Define types
687  using namespace Indices;
688  using DGSpace = TypeTree::Child<LFSV,_0>;
689  using size_type = typename DGSpace::Traits::SizeType;
690 
691  // Get local function space that is identical for all components
692  const auto& dgspace = child(lfsv,_0);
693 
694  // Reference to cell
695  const auto& cell = eg.entity();
696 
697  // Get geometry
698  auto geo = eg.geometry();
699 
700  // loop over quadrature points
701  const int order_s = dgspace.finiteElement().localBasis().order();
702  const int intorder = overintegration+2*order_s;
703  for (const auto &qp : quadratureRule(geo,intorder))
704  {
705  // evaluate right hand side
706  auto j = param.j(cell,qp.position());
707 
708  // evaluate basis functions
709  const auto& phi = cache[order_s].evaluateFunction(qp.position(),dgspace.finiteElement().localBasis());
710 
711  // integrate
712  auto factor = qp.weight() * geo.integrationElement(qp.position());
713  for (size_type k=0; k<dim*2; k++) // for all components
714  for (size_type i=0; i<dgspace.size(); i++) // for all test functions of this component
715  r.accumulate(lfsv.child(k),i,-j[k]*phi[i]*factor);
716  }
717  }
718 
720  void setTime (typename T::Traits::RangeFieldType t)
721  {
722  }
723 
725  void preStep (typename T::Traits::RangeFieldType time, typename T::Traits::RangeFieldType dt,
726  int stages)
727  {
728  }
729 
731  void preStage (typename T::Traits::RangeFieldType time, int r)
732  {
733  }
734 
736  void postStage ()
737  {
738  }
739 
741  typename T::Traits::RangeFieldType suggestTimestep (typename T::Traits::RangeFieldType dt) const
742  {
743  return dt;
744  }
745 
746  private:
747  T& param;
748  int overintegration;
749  typedef typename FEM::Traits::FiniteElementType::Traits::LocalBasisType LocalBasisType;
751  std::vector<Cache> cache;
752  };
753 
754 
755 
767  template<typename T, typename FEM>
769  public NumericalJacobianApplyVolume<DGMaxwellTemporalOperator<T,FEM> >,
771  public InstationaryLocalOperatorDefaultMethods<typename T::Traits::RangeFieldType>
772  {
773  enum { dim = T::Traits::GridViewType::dimension };
774  public:
775  // pattern assembly flags
776  enum { doPatternVolume = true };
777 
778  // residual assembly flags
779  enum { doAlphaVolume = true };
780 
781  DGMaxwellTemporalOperator (T& param_, int overintegration_=0)
782  : param(param_), overintegration(overintegration_), cache(20)
783  {}
784 
785  // define sparsity pattern of operator representation
786  template<typename LFSU, typename LFSV, typename LocalPattern>
787  void pattern_volume (const LFSU& lfsu, const LFSV& lfsv,
788  LocalPattern& pattern) const
789  {
790  // paranoia check number of number of components
792  static_assert(TypeTree::StaticDegree<LFSV>::value==dim*2, "need exactly dim*2 components!");
793 
794  for (size_t k=0; k<TypeTree::degree(lfsv); k++)
795  for (size_t i=0; i<lfsv.child(k).size(); ++i)
796  for (size_t j=0; j<lfsu.child(k).size(); ++j)
797  pattern.addLink(lfsv.child(k),i,lfsu.child(k),j);
798  }
799 
800  // volume integral depending on test and ansatz functions
801  template<typename EG, typename LFSU, typename X, typename LFSV, typename R>
802  void alpha_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv, R& r) const
803  {
804  // Define types
805  using namespace Indices;
806  using DGSpace = TypeTree::Child<LFSV,_0>;
807  using RF = typename DGSpace::Traits::FiniteElementType::
808  Traits::LocalBasisType::Traits::RangeFieldType;
809  using size_type = typename DGSpace::Traits::SizeType;
810 
811  // get local function space that is identical for all components
812  const auto& dgspace = child(lfsv,_0);
813 
814  // Get geometry
815  auto geo = eg.geometry();
816 
817  // Initialize vectors outside for loop
818  Dune::FieldVector<RF,dim*2> u(0.0);
819 
820  // loop over quadrature points
821  const int order = dgspace.finiteElement().localBasis().order();
822  const int intorder = overintegration+2*order;
823  for (const auto& qp : quadratureRule(geo,intorder))
824  {
825  // evaluate basis functions
826  const auto& phi = cache[order].evaluateFunction(qp.position(),dgspace.finiteElement().localBasis());
827 
828  // evaluate u
829  for (size_type k=0; k<dim*2; k++){ // for all components
830  u[k] = 0.0;
831  for (size_type j=0; j<dgspace.size(); j++) // for all basis functions
832  u[k] += x(lfsv.child(k),j)*phi[j];
833  }
834 
835  // integrate
836  auto factor = qp.weight() * geo.integrationElement(qp.position());
837  for (size_type k=0; k<dim*2; k++) // for all components
838  for (size_type i=0; i<dgspace.size(); i++) // for all test functions of this component
839  r.accumulate(lfsv.child(k),i,u[k]*phi[i]*factor);
840  }
841  }
842 
843  // jacobian of volume term
844  template<typename EG, typename LFSU, typename X, typename LFSV, typename M>
845  void jacobian_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv,
846  M& mat) const
847  {
848  // Define types
849  using namespace Indices;
850  using DGSpace = TypeTree::Child<LFSV,_0>;
851  using size_type = typename DGSpace::Traits::SizeType;
852 
853  // Get local function space that is identical for all components
854  const auto& dgspace = child(lfsv,_0);
855 
856  // Get geometry
857  auto geo = eg.geometry();
858 
859  // Loop over quadrature points
860  const int order = dgspace.finiteElement().localBasis().order();
861  const int intorder = overintegration+2*order;
862  for (const auto& qp : quadratureRule(geo,intorder))
863  {
864  // Evaluate basis functions
865  const auto& phi = cache[order].evaluateFunction(qp.position(),dgspace.finiteElement().localBasis());
866 
867  // Integrate
868  auto factor = qp.weight() * geo.integrationElement(qp.position());
869  for (size_type k=0; k<dim*2; k++) // for all components
870  for (size_type i=0; i<dgspace.size(); i++) // for all test functions of this component
871  for (size_type j=0; j<dgspace.size(); j++) // for all ansatz functions of this component
872  mat.accumulate(lfsv.child(k),i,lfsu.child(k),j,phi[j]*phi[i]*factor);
873  }
874  }
875 
876  private:
877  T& param;
878  int overintegration;
879  typedef typename FEM::Traits::FiniteElementType::Traits::LocalBasisType LocalBasisType;
881  std::vector<Cache> cache;
882  };
883 
884  }
885 }
886 
887 #endif // DUNE_PDELAB_LOCALOPERATOR_MAXWELLDG_HH
QuadratureRuleWrapper< QuadratureRule< typename Geometry::ctype, Geometry::mydimension > > quadratureRule(const Geometry &geo, std::size_t order, QuadratureType::Enum quadrature_type=QuadratureType::GaussLegendre)
Returns a quadrature rule for the given geometry.
Definition: quadraturerules.hh:117
T::Traits::RangeFieldType suggestTimestep(typename T::Traits::RangeFieldType dt) const
to be called once before each stage
Definition: maxwelldg.hh:741
DGMaxwellTemporalOperator(T &param_, int overintegration_=0)
Definition: maxwelldg.hh:781
const IG & ig
Definition: constraints.hh:148
static const int dim
Definition: adaptivity.hh:83
void setTime(typename T::Traits::RangeFieldType t)
set time in parameter class
Definition: maxwelldg.hh:720
Default class for additional methods in instationary local operators.
Definition: idefault.hh:89
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
void alpha_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
Definition: maxwelldg.hh:802
sparsity pattern generator
Definition: pattern.hh:29
const std::string s
Definition: function.hh:1101
Implement jacobian_skeleton() based on alpha_skeleton()
Definition: numericaljacobian.hh:156
Implements linear and nonlinear versions of jacobian_apply_boundary() based on alpha_boundary() ...
Definition: numericaljacobianapply.hh:285
Implements linear and nonlinear versions of jacobian_apply_volume() based on alpha_volume() ...
Definition: numericaljacobianapply.hh:32
void lambda_volume(const EG &eg, const LFSV &lfsv, R &r) const
Definition: maxwelldg.hh:684
void pattern_volume(const LFSU &lfsu, const LFSV &lfsv, LocalPattern &pattern) const
Definition: maxwelldg.hh:787
void postStage()
to be called once at the end of each stage
Definition: maxwelldg.hh:736
DGMaxwellSpatialOperator(T &param_, int overintegration_=0)
Definition: maxwelldg.hh:330
Implements linear and nonlinear versions of jacobian_apply_skeleton() based on alpha_skeleton() ...
Definition: numericaljacobianapply.hh:180
store values of basis functions and gradients in a cache
Definition: localbasiscache.hh:17
void jacobian_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, M &mat) const
Definition: maxwelldg.hh:845
static void eigenvectors(T1 eps, T1 mu, const Dune::FieldVector< T2, dim > &n, Dune::FieldMatrix< T3, 2 *dim, 2 *dim > &R)
Definition: maxwelldg.hh:76
void preStep(typename T::Traits::RangeFieldType time, typename T::Traits::RangeFieldType dt, int stages)
to be called once before each time step
Definition: maxwelldg.hh:725
Definition: maxwelldg.hh:768
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
void alpha_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
Definition: maxwelldg.hh:337
sparsity pattern generator
Definition: pattern.hh:13
void preStage(typename T::Traits::RangeFieldType time, int r)
to be called once before each stage
Definition: maxwelldg.hh:731
const Entity & e
Definition: localfunctionspace.hh:111
static void eigenvalues(T1 eps, T1 mu, const Dune::FieldVector< T2, 2 *dim > &e)
Definition: maxwelldg.hh:53
void alpha_boundary(const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, R &r_s) const
Definition: maxwelldg.hh:574
Implement jacobian_boundary() based on alpha_boundary()
Definition: numericaljacobian.hh:250
Definition: maxwelldg.hh:29
Implement jacobian_volume() based on alpha_volume()
Definition: numericaljacobian.hh:31
void alpha_skeleton(const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const LFSV &lfsv_n, R &r_s, R &r_n) const
Definition: maxwelldg.hh:437
ReferenceElementWrapper< ReferenceElement< typename Geometry::ctype, Geometry::mydimension > > referenceElement(const Geometry &geo)
Returns the reference element for the given geometry.
Definition: referenceelements.hh:144
Default flags for all local operators.
Definition: flags.hh:18
Definition: maxwelldg.hh:304