MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Implementation of Bilinear Form Integrators
13 
14 #include "fem.hpp"
15 #include <cmath>
16 #include <algorithm>
17 
18 using namespace std;
19 
20 namespace mfem
21 {
22 
23 void BilinearFormIntegrator::AssembleElementMatrix (
25  DenseMatrix &elmat )
26 {
27  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix (...)\n"
28  " is not implemented fot this class.");
29 }
30 
31 void BilinearFormIntegrator::AssembleElementMatrix2 (
32  const FiniteElement &el1, const FiniteElement &el2,
34 {
35  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix2 (...)\n"
36  " is not implemented fot this class.");
37 }
38 
39 void BilinearFormIntegrator::AssembleFaceMatrix (
40  const FiniteElement &el1, const FiniteElement &el2,
42 {
43  mfem_error ("BilinearFormIntegrator::AssembleFaceMatrix (...)\n"
44  " is not implemented fot this class.");
45 }
46 
47 void BilinearFormIntegrator::AssembleFaceMatrix(
48  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
50  DenseMatrix &elmat)
51 {
52  MFEM_ABORT("AssembleFaceMatrix (mixed form) is not implemented for this"
53  " Integrator class.");
54 }
55 
56 void BilinearFormIntegrator::AssembleElementVector(
57  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
58  Vector &elvect)
59 {
60  mfem_error("BilinearFormIntegrator::AssembleElementVector\n"
61  " is not implemented fot this class.");
62 }
63 
64 
65 void TransposeIntegrator::AssembleElementMatrix (
67 {
68  bfi -> AssembleElementMatrix (el, Trans, bfi_elmat);
69  // elmat = bfi_elmat^t
70  elmat.Transpose (bfi_elmat);
71 }
72 
73 void TransposeIntegrator::AssembleElementMatrix2 (
74  const FiniteElement &trial_fe, const FiniteElement &test_fe,
76 {
77  bfi -> AssembleElementMatrix2 (test_fe, trial_fe, Trans, bfi_elmat);
78  // elmat = bfi_elmat^t
79  elmat.Transpose (bfi_elmat);
80 }
81 
82 void TransposeIntegrator::AssembleFaceMatrix (
83  const FiniteElement &el1, const FiniteElement &el2,
85 {
86  bfi -> AssembleFaceMatrix (el1, el2, Trans, bfi_elmat);
87  // elmat = bfi_elmat^t
88  elmat.Transpose (bfi_elmat);
89 }
90 
91 void LumpedIntegrator::AssembleElementMatrix (
93 {
94  bfi -> AssembleElementMatrix (el, Trans, elmat);
95  elmat.Lump();
96 }
97 
98 void InverseIntegrator::AssembleElementMatrix(
100 {
101  integrator->AssembleElementMatrix(el, Trans, elmat);
102  elmat.Invert();
103 }
104 
105 void SumIntegrator::AssembleElementMatrix(
107 {
108  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
109 
110  integrators[0]->AssembleElementMatrix(el, Trans, elmat);
111  for (int i = 1; i < integrators.Size(); i++)
112  {
113  integrators[i]->AssembleElementMatrix(el, Trans, elem_mat);
114  elmat += elem_mat;
115  }
116 }
117 
118 SumIntegrator::~SumIntegrator()
119 {
120  if (own_integrators)
121  {
122  for (int i = 0; i < integrators.Size(); i++)
123  {
124  delete integrators[i];
125  }
126  }
127 }
128 
129 void MixedScalarIntegrator::AssembleElementMatrix2(
130  const FiniteElement &trial_fe, const FiniteElement &test_fe,
132 {
133  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
134  this->FiniteElementTypeFailureMessage());
135 
136  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
137  bool same_shapes = same_calc_shape && (&trial_fe == &test_fe);
138 
139 #ifdef MFEM_THREAD_SAFE
140  Vector test_shape(test_nd);
141  Vector trial_shape;
142 #else
143  test_shape.SetSize(test_nd);
144 #endif
145  if (same_shapes)
146  {
147  trial_shape.NewDataAndSize(test_shape.GetData(), trial_nd);
148  }
149  else
150  {
151  trial_shape.SetSize(trial_nd);
152  }
153 
154  elmat.SetSize(test_nd, trial_nd);
155 
156  const IntegrationRule *ir = IntRule;
157  if (ir == NULL)
158  {
159  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
160  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
161  }
162 
163  elmat = 0.0;
164  for (i = 0; i < ir->GetNPoints(); i++)
165  {
166  const IntegrationPoint &ip = ir->IntPoint(i);
167  Trans.SetIntPoint(&ip);
168 
169  this->CalcTestShape(test_fe, Trans, test_shape);
170  this->CalcTrialShape(trial_fe, Trans, trial_shape);
171 
172  double w = Trans.Weight() * ip.weight;
173 
174  if (Q)
175  {
176  w *= Q->Eval(Trans, ip);
177  }
178  AddMult_a_VWt(w, test_shape, trial_shape, elmat);
179  }
180 #ifndef MFEM_THREAD_SAFE
181  if (same_shapes)
182  {
183  trial_shape.SetDataAndSize(NULL, 0);
184  }
185 #endif
186 }
187 
188 void MixedVectorIntegrator::AssembleElementMatrix2(
189  const FiniteElement &trial_fe, const FiniteElement &test_fe,
191 {
192  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
193  this->FiniteElementTypeFailureMessage());
194 
195  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
196  int spaceDim = Trans.GetSpaceDim();
197  bool same_shapes = same_calc_shape && (&trial_fe == &test_fe);
198 
199 #ifdef MFEM_THREAD_SAFE
200  Vector V(VQ ? VQ->GetVDim() : 0);
201  Vector D(DQ ? DQ->GetVDim() : 0);
202  DenseMatrix M(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
203  DenseMatrix test_shape(test_nd, spaceDim);
204  DenseMatrix trial_shape;
205  DenseMatrix test_shape_tmp(test_nd, spaceDim);
206 #else
207  V.SetSize(VQ ? VQ->GetVDim() : 0);
208  D.SetSize(DQ ? DQ->GetVDim() : 0);
209  M.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
210  test_shape.SetSize(test_nd, spaceDim);
211  test_shape_tmp.SetSize(test_nd, spaceDim);
212 #endif
213  if (same_shapes)
214  {
215  trial_shape.Reset(test_shape.Data(), trial_nd, spaceDim);
216  }
217  else
218  {
219  trial_shape.SetSize(trial_nd, spaceDim);
220  }
221 
222  elmat.SetSize(test_nd, trial_nd);
223 
224  const IntegrationRule *ir = IntRule;
225  if (ir == NULL)
226  {
227  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
228  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
229  }
230 
231  elmat = 0.0;
232  for (i = 0; i < ir->GetNPoints(); i++)
233  {
234  const IntegrationPoint &ip = ir->IntPoint(i);
235  Trans.SetIntPoint(&ip);
236 
237  this->CalcTestShape(test_fe, Trans, test_shape);
238  if (!same_shapes)
239  {
240  this->CalcTrialShape(trial_fe, Trans, trial_shape);
241  }
242 
243  double w = Trans.Weight() * ip.weight;
244 
245  if (MQ)
246  {
247  MQ->Eval(M, Trans, ip);
248  M *= w;
249  Mult(test_shape, M, test_shape_tmp);
250  AddMultABt(test_shape_tmp, trial_shape, elmat);
251  }
252  else if (DQ)
253  {
254  DQ->Eval(D, Trans, ip);
255  D *= w;
256  AddMultADBt(test_shape, D, trial_shape, elmat);
257  }
258  else if (VQ)
259  {
260  VQ->Eval(V, Trans, ip);
261  V *= w;
262  for (int j=0; j<test_nd; j++)
263  {
264  test_shape_tmp(j,0) = test_shape(j,1) * V(2) -
265  test_shape(j,2) * V(1);
266  test_shape_tmp(j,1) = test_shape(j,2) * V(0) -
267  test_shape(j,0) * V(2);
268  test_shape_tmp(j,2) = test_shape(j,0) * V(1) -
269  test_shape(j,1) * V(0);
270  }
271  AddMultABt(test_shape_tmp, trial_shape, elmat);
272  }
273  else
274  {
275  if (Q)
276  {
277  w *= Q -> Eval (Trans, ip);
278  }
279  if (same_shapes)
280  {
281  AddMult_a_AAt (w, test_shape, elmat);
282  }
283  else
284  {
285  AddMult_a_ABt (w, test_shape, trial_shape, elmat);
286  }
287  }
288  }
289 #ifndef MFEM_THREAD_SAFE
290  if (same_shapes)
291  {
292  trial_shape.ClearExternalData();
293  }
294 #endif
295 }
296 
297 void MixedScalarVectorIntegrator::AssembleElementMatrix2(
298  const FiniteElement &trial_fe, const FiniteElement &test_fe,
300 {
301  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
302  this->FiniteElementTypeFailureMessage());
303 
304  const FiniteElement * vec_fe = transpose?&trial_fe:&test_fe;
305  const FiniteElement * sca_fe = transpose?&test_fe:&trial_fe;
306 
307  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
308  int sca_nd = sca_fe->GetDof();
309  int vec_nd = vec_fe->GetDof();
310  int spaceDim = Trans.GetSpaceDim();
311  double vtmp;
312 
313 #ifdef MFEM_THREAD_SAFE
314  Vector V(VQ ? VQ->GetVDim() : 0);
315  DenseMatrix vshape(vec_nd, spaceDim);
316  Vector shape(sca_nd);
317  Vector vshape_tmp(vec_nd);
318 #else
319  V.SetSize(VQ ? VQ->GetVDim() : 0);
320  vshape.SetSize(vec_nd, spaceDim);
321  shape.SetSize(sca_nd);
322  vshape_tmp.SetSize(vec_nd);
323 #endif
324 
325  Vector V_test(transpose?shape.GetData():vshape_tmp.GetData(),test_nd);
326  Vector W_trial(transpose?vshape_tmp.GetData():shape.GetData(),trial_nd);
327 
328  elmat.SetSize(test_nd, trial_nd);
329 
330  const IntegrationRule *ir = IntRule;
331  if (ir == NULL)
332  {
333  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
334  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
335  }
336 
337  elmat = 0.0;
338  for (i = 0; i < ir->GetNPoints(); i++)
339  {
340  const IntegrationPoint &ip = ir->IntPoint(i);
341  Trans.SetIntPoint(&ip);
342 
343  this->CalcShape(*sca_fe, Trans, shape);
344  this->CalcVShape(*vec_fe, Trans, vshape);
345 
346  double w = Trans.Weight() * ip.weight;
347 
348  VQ->Eval(V, Trans, ip);
349  V *= w;
350 
351  if ( vec_fe->GetDim() == 2 && cross_2d )
352  {
353  vtmp = V[0];
354  V[0] = -V[1];
355  V[1] = vtmp;
356  }
357 
358  vshape.Mult(V,vshape_tmp);
359 
360  AddMultVWt(V_test, W_trial, elmat);
361  }
362 }
363 
364 void DiffusionIntegrator::AssembleElementMatrix
366  DenseMatrix &elmat )
367 {
368  int nd = el.GetDof();
369  int dim = el.GetDim();
370  int spaceDim = Trans.GetSpaceDim();
371  bool square = (dim == spaceDim);
372  double w;
373 
374 #ifdef MFEM_THREAD_SAFE
375  DenseMatrix dshape(nd,dim), dshapedxt(nd,spaceDim), invdfdx(dim,spaceDim);
376 #else
377  dshape.SetSize(nd,dim);
378  dshapedxt.SetSize(nd,spaceDim);
379  invdfdx.SetSize(dim,spaceDim);
380 #endif
381  elmat.SetSize(nd);
382 
383  const IntegrationRule *ir = IntRule;
384  if (ir == NULL)
385  {
386  int order;
387  if (el.Space() == FunctionSpace::Pk)
388  {
389  order = 2*el.GetOrder() - 2;
390  }
391  else
392  // order = 2*el.GetOrder() - 2; // <-- this seems to work fine too
393  {
394  order = 2*el.GetOrder() + dim - 1;
395  }
396 
397  if (el.Space() == FunctionSpace::rQk)
398  {
399  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
400  }
401  else
402  {
403  ir = &IntRules.Get(el.GetGeomType(), order);
404  }
405  }
406 
407  elmat = 0.0;
408  for (int i = 0; i < ir->GetNPoints(); i++)
409  {
410  const IntegrationPoint &ip = ir->IntPoint(i);
411  el.CalcDShape(ip, dshape);
412 
413  Trans.SetIntPoint(&ip);
414  w = Trans.Weight();
415  w = ip.weight / (square ? w : w*w*w);
416  // AdjugateJacobian = / adj(J), if J is square
417  // \ adj(J^t.J).J^t, otherwise
418  Mult(dshape, Trans.AdjugateJacobian(), dshapedxt);
419  if (!MQ)
420  {
421  if (Q)
422  {
423  w *= Q->Eval(Trans, ip);
424  }
425  AddMult_a_AAt(w, dshapedxt, elmat);
426  }
427  else
428  {
429  MQ->Eval(invdfdx, Trans, ip);
430  invdfdx *= w;
431  Mult(dshapedxt, invdfdx, dshape);
432  AddMultABt(dshape, dshapedxt, elmat);
433  }
434  }
435 }
436 
437 void DiffusionIntegrator::AssembleElementMatrix2(
438  const FiniteElement &trial_fe, const FiniteElement &test_fe,
439  ElementTransformation &Trans, DenseMatrix &elmat)
440 {
441  int tr_nd = trial_fe.GetDof();
442  int te_nd = test_fe.GetDof();
443  int dim = trial_fe.GetDim();
444  int spaceDim = Trans.GetSpaceDim();
445  bool square = (dim == spaceDim);
446  double w;
447 
448 #ifdef MFEM_THREAD_SAFE
449  DenseMatrix dshape(tr_nd, dim), dshapedxt(tr_nd, spaceDim);
450  DenseMatrix te_dshape(te_nd, dim), te_dshapedxt(te_nd, spaceDim);
451  DenseMatrix invdfdx(dim, spaceDim);
452 #else
453  dshape.SetSize(tr_nd, dim);
454  dshapedxt.SetSize(tr_nd, spaceDim);
455  te_dshape.SetSize(te_nd, dim);
456  te_dshapedxt.SetSize(te_nd, spaceDim);
457  invdfdx.SetSize(dim, spaceDim);
458 #endif
459  elmat.SetSize(te_nd, tr_nd);
460 
461  const IntegrationRule *ir = IntRule;
462  if (ir == NULL)
463  {
464  int order;
465  if (trial_fe.Space() == FunctionSpace::Pk)
466  {
467  order = trial_fe.GetOrder() + test_fe.GetOrder() - 2;
468  }
469  else
470  {
471  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim - 1;
472  }
473 
474  if (trial_fe.Space() == FunctionSpace::rQk)
475  {
476  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
477  }
478  else
479  {
480  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
481  }
482  }
483 
484  elmat = 0.0;
485  for (int i = 0; i < ir->GetNPoints(); i++)
486  {
487  const IntegrationPoint &ip = ir->IntPoint(i);
488  trial_fe.CalcDShape(ip, dshape);
489  test_fe.CalcDShape(ip, te_dshape);
490 
491  Trans.SetIntPoint(&ip);
492  CalcAdjugate(Trans.Jacobian(), invdfdx);
493  w = Trans.Weight();
494  w = ip.weight / (square ? w : w*w*w);
495  Mult(dshape, invdfdx, dshapedxt);
496  Mult(te_dshape, invdfdx, te_dshapedxt);
497  // invdfdx, dshape, and te_dshape no longer needed
498  if (!MQ)
499  {
500  if (Q)
501  {
502  w *= Q->Eval(Trans, ip);
503  }
504  dshapedxt *= w;
505  AddMultABt(te_dshapedxt, dshapedxt, elmat);
506  }
507  else
508  {
509  MQ->Eval(invdfdx, Trans, ip);
510  invdfdx *= w;
511  Mult(te_dshapedxt, invdfdx, te_dshape);
512  AddMultABt(te_dshape, dshapedxt, elmat);
513  }
514  }
515 }
516 
517 void DiffusionIntegrator::AssembleElementVector(
518  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
519  Vector &elvect)
520 {
521  int nd = el.GetDof();
522  int dim = el.GetDim();
523  double w;
524 
525 #ifdef MFEM_THREAD_SAFE
526  DenseMatrix dshape(nd,dim), invdfdx(dim), mq(dim);
527 #else
528  dshape.SetSize(nd,dim);
529  invdfdx.SetSize(dim);
530  mq.SetSize(dim);
531 #endif
532  vec.SetSize(dim);
533  pointflux.SetSize(dim);
534 
535  elvect.SetSize(nd);
536 
537  const IntegrationRule *ir = IntRule;
538  if (ir == NULL)
539  {
540  int order;
541  if (el.Space() == FunctionSpace::Pk)
542  {
543  order = 2*el.GetOrder() - 2;
544  }
545  else
546  // order = 2*el.GetOrder() - 2; // <-- this seems to work fine too
547  {
548  order = 2*el.GetOrder() + dim - 1;
549  }
550 
551  if (el.Space() == FunctionSpace::rQk)
552  {
553  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
554  }
555  else
556  {
557  ir = &IntRules.Get(el.GetGeomType(), order);
558  }
559  }
560 
561  elvect = 0.0;
562  for (int i = 0; i < ir->GetNPoints(); i++)
563  {
564  const IntegrationPoint &ip = ir->IntPoint(i);
565  el.CalcDShape(ip, dshape);
566 
567  Tr.SetIntPoint(&ip);
568  CalcAdjugate(Tr.Jacobian(), invdfdx); // invdfdx = adj(J)
569  w = ip.weight / Tr.Weight();
570 
571  if (!MQ)
572  {
573  dshape.MultTranspose(elfun, vec);
574  invdfdx.MultTranspose(vec, pointflux);
575  if (Q)
576  {
577  w *= Q->Eval(Tr, ip);
578  }
579  }
580  else
581  {
582 
583  dshape.MultTranspose(elfun, pointflux);
584  invdfdx.MultTranspose(pointflux, vec);
585  MQ->Eval(mq, Tr, ip);
586  mq.Mult(vec, pointflux);
587  }
588  pointflux *= w;
589  invdfdx.Mult(pointflux, vec);
590  dshape.AddMult(vec, elvect);
591  }
592 }
593 
594 void DiffusionIntegrator::ComputeElementFlux
596  Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef )
597 {
598  int i, j, nd, dim, spaceDim, fnd;
599 
600  nd = el.GetDof();
601  dim = el.GetDim();
602  spaceDim = Trans.GetSpaceDim();
603 
604 #ifdef MFEM_THREAD_SAFE
605  DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim);
606 #else
607  dshape.SetSize(nd,dim);
608  invdfdx.SetSize(dim, spaceDim);
609 #endif
610  vec.SetSize(dim);
611  pointflux.SetSize(spaceDim);
612 
613  const IntegrationRule &ir = fluxelem.GetNodes();
614  fnd = ir.GetNPoints();
615  flux.SetSize( fnd * spaceDim );
616 
617  for (i = 0; i < fnd; i++)
618  {
619  const IntegrationPoint &ip = ir.IntPoint(i);
620  el.CalcDShape(ip, dshape);
621  dshape.MultTranspose(u, vec);
622 
623  Trans.SetIntPoint (&ip);
624  CalcInverse(Trans.Jacobian(), invdfdx);
625  invdfdx.MultTranspose(vec, pointflux);
626 
627  if (!MQ)
628  {
629  if (Q && with_coef)
630  {
631  pointflux *= Q->Eval(Trans,ip);
632  }
633  for (j = 0; j < spaceDim; j++)
634  {
635  flux(fnd*j+i) = pointflux(j);
636  }
637  }
638  else
639  {
640  // assuming dim == spaceDim
641  MFEM_ASSERT(dim == spaceDim, "TODO");
642  MQ->Eval(invdfdx, Trans, ip);
643  invdfdx.Mult(pointflux, vec);
644  for (j = 0; j < dim; j++)
645  {
646  flux(fnd*j+i) = vec(j);
647  }
648  }
649  }
650 }
651 
652 double DiffusionIntegrator::ComputeFluxEnergy
653 ( const FiniteElement &fluxelem, ElementTransformation &Trans,
654  Vector &flux, Vector* d_energy)
655 {
656  int nd = fluxelem.GetDof();
657  int dim = fluxelem.GetDim();
658  int spaceDim = Trans.GetSpaceDim();
659 
660 #ifdef MFEM_THREAD_SAFE
661  DenseMatrix mq;
662 #endif
663 
664  shape.SetSize(nd);
665  pointflux.SetSize(spaceDim);
666  if (d_energy) { vec.SetSize(dim); }
667  if (MQ) { mq.SetSize(dim); }
668 
669  int order = 2 * fluxelem.GetOrder(); // <--
670  const IntegrationRule *ir = &IntRules.Get(fluxelem.GetGeomType(), order);
671 
672  double energy = 0.0;
673  if (d_energy) { *d_energy = 0.0; }
674 
675  for (int i = 0; i < ir->GetNPoints(); i++)
676  {
677  const IntegrationPoint &ip = ir->IntPoint(i);
678  fluxelem.CalcShape(ip, shape);
679 
680  pointflux = 0.0;
681  for (int k = 0; k < spaceDim; k++)
682  {
683  for (int j = 0; j < nd; j++)
684  {
685  pointflux(k) += flux(k*nd+j)*shape(j);
686  }
687  }
688 
689  Trans.SetIntPoint(&ip);
690  double w = Trans.Weight() * ip.weight;
691 
692  if (!MQ)
693  {
694  double e = (pointflux * pointflux);
695  if (Q) { e *= Q->Eval(Trans, ip); }
696  energy += w * e;
697  }
698  else
699  {
700  MQ->Eval(mq, Trans, ip);
701  energy += w * mq.InnerProduct(pointflux, pointflux);
702  }
703 
704  if (d_energy)
705  {
706  // transform pointflux to the ref. domain and integrate the components
707  Trans.Jacobian().MultTranspose(pointflux, vec);
708  for (int k = 0; k < dim; k++)
709  {
710  (*d_energy)[k] += w * vec[k] * vec[k];
711  }
712  // TODO: Q, MQ
713  }
714  }
715 
716  return energy;
717 }
718 
719 
720 void MassIntegrator::AssembleElementMatrix
722  DenseMatrix &elmat )
723 {
724  int nd = el.GetDof();
725  // int dim = el.GetDim();
726  double w;
727 
728  elmat.SetSize(nd);
729  shape.SetSize(nd);
730 
731  const IntegrationRule *ir = IntRule;
732  if (ir == NULL)
733  {
734  // int order = 2 * el.GetOrder();
735  int order = 2 * el.GetOrder() + Trans.OrderW();
736 
737  if (el.Space() == FunctionSpace::rQk)
738  {
739  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
740  }
741  else
742  {
743  ir = &IntRules.Get(el.GetGeomType(), order);
744  }
745  }
746 
747  elmat = 0.0;
748  for (int i = 0; i < ir->GetNPoints(); i++)
749  {
750  const IntegrationPoint &ip = ir->IntPoint(i);
751  el.CalcShape(ip, shape);
752 
753  Trans.SetIntPoint (&ip);
754  w = Trans.Weight() * ip.weight;
755  if (Q)
756  {
757  w *= Q -> Eval(Trans, ip);
758  }
759 
760  AddMult_a_VVt(w, shape, elmat);
761  }
762 }
763 
764 void MassIntegrator::AssembleElementMatrix2(
765  const FiniteElement &trial_fe, const FiniteElement &test_fe,
766  ElementTransformation &Trans, DenseMatrix &elmat)
767 {
768  int tr_nd = trial_fe.GetDof();
769  int te_nd = test_fe.GetDof();
770  // int dim = trial_fe.GetDim();
771  double w;
772 
773  elmat.SetSize (te_nd, tr_nd);
774  shape.SetSize (tr_nd);
775  te_shape.SetSize (te_nd);
776 
777  const IntegrationRule *ir = IntRule;
778  if (ir == NULL)
779  {
780  int order = trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW();
781 
782  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
783  }
784 
785  elmat = 0.0;
786  for (int i = 0; i < ir->GetNPoints(); i++)
787  {
788  const IntegrationPoint &ip = ir->IntPoint(i);
789  trial_fe.CalcShape(ip, shape);
790  test_fe.CalcShape(ip, te_shape);
791 
792  Trans.SetIntPoint (&ip);
793  w = Trans.Weight() * ip.weight;
794  if (Q)
795  {
796  w *= Q -> Eval(Trans, ip);
797  }
798 
799  te_shape *= w;
800  AddMultVWt(te_shape, shape, elmat);
801  }
802 }
803 
804 
805 void BoundaryMassIntegrator::AssembleFaceMatrix(
806  const FiniteElement &el1, const FiniteElement &el2,
808 {
809  MFEM_ASSERT(Trans.Elem2No < 0,
810  "support for interior faces is not implemented");
811 
812  int nd1 = el1.GetDof();
813  double w;
814 
815  elmat.SetSize(nd1);
816  shape.SetSize(nd1);
817 
818  const IntegrationRule *ir = IntRule;
819  if (ir == NULL)
820  {
821  int order = 2 * el1.GetOrder();
822 
823  ir = &IntRules.Get(Trans.FaceGeom, order);
824  }
825 
826  elmat = 0.0;
827  for (int i = 0; i < ir->GetNPoints(); i++)
828  {
829  const IntegrationPoint &ip = ir->IntPoint(i);
830  IntegrationPoint eip;
831  Trans.Loc1.Transform(ip, eip);
832  el1.CalcShape(eip, shape);
833 
834  Trans.Face->SetIntPoint(&ip);
835  w = Trans.Face->Weight() * ip.weight;
836  if (Q)
837  {
838  w *= Q -> Eval(*Trans.Face, ip);
839  }
840 
841  AddMult_a_VVt(w, shape, elmat);
842  }
843 }
844 
845 
846 void ConvectionIntegrator::AssembleElementMatrix(
847  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
848 {
849  int nd = el.GetDof();
850  int dim = el.GetDim();
851 
852  elmat.SetSize(nd);
853  dshape.SetSize(nd,dim);
854  adjJ.SetSize(dim);
855  shape.SetSize(nd);
856  vec2.SetSize(dim);
857  BdFidxT.SetSize(nd);
858 
859  Vector vec1;
860 
861  const IntegrationRule *ir = IntRule;
862  if (ir == NULL)
863  {
864  int order = Trans.OrderGrad(&el) + Trans.Order() + el.GetOrder();
865  ir = &IntRules.Get(el.GetGeomType(), order);
866  }
867 
868  Q.Eval(Q_ir, Trans, *ir);
869 
870  elmat = 0.0;
871  for (int i = 0; i < ir->GetNPoints(); i++)
872  {
873  const IntegrationPoint &ip = ir->IntPoint(i);
874  el.CalcDShape(ip, dshape);
875  el.CalcShape(ip, shape);
876 
877  Trans.SetIntPoint(&ip);
878  CalcAdjugate(Trans.Jacobian(), adjJ);
879  Q_ir.GetColumnReference(i, vec1);
880  vec1 *= alpha * ip.weight;
881 
882  adjJ.Mult(vec1, vec2);
883  dshape.Mult(vec2, BdFidxT);
884 
885  AddMultVWt(shape, BdFidxT, elmat);
886  }
887 }
888 
889 
890 void GroupConvectionIntegrator::AssembleElementMatrix(
891  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
892 {
893  int nd = el.GetDof();
894  int dim = el.GetDim();
895 
896  elmat.SetSize(nd);
897  dshape.SetSize(nd,dim);
898  adjJ.SetSize(dim);
899  shape.SetSize(nd);
900  grad.SetSize(nd,dim);
901 
902  const IntegrationRule *ir = IntRule;
903  if (ir == NULL)
904  {
905  int order = Trans.OrderGrad(&el) + el.GetOrder();
906  ir = &IntRules.Get(el.GetGeomType(), order);
907  }
908 
909  Q.Eval(Q_nodal, Trans, el.GetNodes()); // sets the size of Q_nodal
910 
911  elmat = 0.0;
912  for (int i = 0; i < ir->GetNPoints(); i++)
913  {
914  const IntegrationPoint &ip = ir->IntPoint(i);
915  el.CalcDShape(ip, dshape);
916  el.CalcShape(ip, shape);
917 
918  Trans.SetIntPoint(&ip);
919  CalcAdjugate(Trans.Jacobian(), adjJ);
920 
921  Mult(dshape, adjJ, grad);
922 
923  double w = alpha * ip.weight;
924 
925  // elmat(k,l) += \sum_s w*shape(k)*Q_nodal(s,k)*grad(l,s)
926  for (int k = 0; k < nd; k++)
927  {
928  double wsk = w*shape(k);
929  for (int l = 0; l < nd; l++)
930  {
931  double a = 0.0;
932  for (int s = 0; s < dim; s++)
933  {
934  a += Q_nodal(s,k)*grad(l,s);
935  }
936  elmat(k,l) += wsk*a;
937  }
938  }
939  }
940 }
941 
942 
943 void VectorMassIntegrator::AssembleElementMatrix
945  DenseMatrix &elmat )
946 {
947  int nd = el.GetDof();
948  int spaceDim = Trans.GetSpaceDim();
949 
950  double norm;
951 
952  // Get vdim from VQ, MQ, or the space dimension
953  int vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : spaceDim);
954 
955  elmat.SetSize(nd*vdim);
956  shape.SetSize(nd);
957  partelmat.SetSize(nd);
958  if (VQ)
959  {
960  vec.SetSize(vdim);
961  }
962  else if (MQ)
963  {
964  mcoeff.SetSize(vdim);
965  }
966 
967  const IntegrationRule *ir = IntRule;
968  if (ir == NULL)
969  {
970  int order = 2 * el.GetOrder() + Trans.OrderW() + Q_order;
971 
972  if (el.Space() == FunctionSpace::rQk)
973  {
974  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
975  }
976  else
977  {
978  ir = &IntRules.Get(el.GetGeomType(), order);
979  }
980  }
981 
982  elmat = 0.0;
983  for (int s = 0; s < ir->GetNPoints(); s++)
984  {
985  const IntegrationPoint &ip = ir->IntPoint(s);
986  el.CalcShape(ip, shape);
987 
988  Trans.SetIntPoint (&ip);
989  norm = ip.weight * Trans.Weight();
990 
991  MultVVt(shape, partelmat);
992 
993  if (VQ)
994  {
995  VQ->Eval(vec, Trans, ip);
996  for (int k = 0; k < vdim; k++)
997  {
998  elmat.AddMatrix(norm*vec(k), partelmat, nd*k, nd*k);
999  }
1000  }
1001  else if (MQ)
1002  {
1003  MQ->Eval(mcoeff, Trans, ip);
1004  for (int i = 0; i < vdim; i++)
1005  for (int j = 0; j < vdim; j++)
1006  {
1007  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, nd*i, nd*j);
1008  }
1009  }
1010  else
1011  {
1012  if (Q)
1013  {
1014  norm *= Q->Eval(Trans, ip);
1015  }
1016  partelmat *= norm;
1017  for (int k = 0; k < vdim; k++)
1018  {
1019  elmat.AddMatrix(partelmat, nd*k, nd*k);
1020  }
1021  }
1022  }
1023 }
1024 
1025 void VectorMassIntegrator::AssembleElementMatrix2(
1026  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1027  ElementTransformation &Trans, DenseMatrix &elmat)
1028 {
1029  int tr_nd = trial_fe.GetDof();
1030  int te_nd = test_fe.GetDof();
1031  int dim = trial_fe.GetDim();
1032  int vdim;
1033 
1034  double norm;
1035 
1036  // Get vdim from the ElementTransformation Trans ?
1037  vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : (dim));
1038 
1039  elmat.SetSize(te_nd*vdim, tr_nd*vdim);
1040  shape.SetSize(tr_nd);
1041  te_shape.SetSize(te_nd);
1042  partelmat.SetSize(te_nd, tr_nd);
1043  if (VQ)
1044  {
1045  vec.SetSize(vdim);
1046  }
1047  else if (MQ)
1048  {
1049  mcoeff.SetSize(vdim);
1050  }
1051 
1052  const IntegrationRule *ir = IntRule;
1053  if (ir == NULL)
1054  {
1055  int order = (trial_fe.GetOrder() + test_fe.GetOrder() +
1056  Trans.OrderW() + Q_order);
1057 
1058  if (trial_fe.Space() == FunctionSpace::rQk)
1059  {
1060  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1061  }
1062  else
1063  {
1064  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1065  }
1066  }
1067 
1068  elmat = 0.0;
1069  for (int s = 0; s < ir->GetNPoints(); s++)
1070  {
1071  const IntegrationPoint &ip = ir->IntPoint(s);
1072  trial_fe.CalcShape(ip, shape);
1073  test_fe.CalcShape(ip, te_shape);
1074 
1075  Trans.SetIntPoint(&ip);
1076  norm = ip.weight * Trans.Weight();
1077 
1078  MultVWt(te_shape, shape, partelmat);
1079 
1080  if (VQ)
1081  {
1082  VQ->Eval(vec, Trans, ip);
1083  for (int k = 0; k < vdim; k++)
1084  {
1085  elmat.AddMatrix(norm*vec(k), partelmat, te_nd*k, tr_nd*k);
1086  }
1087  }
1088  else if (MQ)
1089  {
1090  MQ->Eval(mcoeff, Trans, ip);
1091  for (int i = 0; i < vdim; i++)
1092  for (int j = 0; j < vdim; j++)
1093  {
1094  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, te_nd*i, tr_nd*j);
1095  }
1096  }
1097  else
1098  {
1099  if (Q)
1100  {
1101  norm *= Q->Eval(Trans, ip);
1102  }
1103  partelmat *= norm;
1104  for (int k = 0; k < vdim; k++)
1105  {
1106  elmat.AddMatrix(partelmat, te_nd*k, tr_nd*k);
1107  }
1108  }
1109  }
1110 }
1111 
1112 void VectorFEDivergenceIntegrator::AssembleElementMatrix2(
1113  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1114  ElementTransformation &Trans, DenseMatrix &elmat)
1115 {
1116  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1117 
1118 #ifdef MFEM_THREAD_SAFE
1119  Vector divshape(trial_nd), shape(test_nd);
1120 #else
1121  divshape.SetSize(trial_nd);
1122  shape.SetSize(test_nd);
1123 #endif
1124 
1125  elmat.SetSize(test_nd, trial_nd);
1126 
1127  const IntegrationRule *ir = IntRule;
1128  if (ir == NULL)
1129  {
1130  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
1131  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1132  }
1133 
1134  elmat = 0.0;
1135  for (i = 0; i < ir->GetNPoints(); i++)
1136  {
1137  const IntegrationPoint &ip = ir->IntPoint(i);
1138  trial_fe.CalcDivShape(ip, divshape);
1139  test_fe.CalcShape(ip, shape);
1140  double w = ip.weight;
1141  if (Q)
1142  {
1143  Trans.SetIntPoint(&ip);
1144  w *= Q->Eval(Trans, ip);
1145  }
1146  shape *= w;
1147  AddMultVWt(shape, divshape, elmat);
1148  }
1149 }
1150 
1151 void VectorFEWeakDivergenceIntegrator::AssembleElementMatrix2(
1152  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1153  ElementTransformation &Trans, DenseMatrix &elmat)
1154 {
1155  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1156  int dim = trial_fe.GetDim();
1157 
1158  MFEM_ASSERT(test_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1159  test_fe.GetMapType() == mfem::FiniteElement::VALUE &&
1161  "Trial space must be H(Curl) and test space must be H_1");
1162 
1163 #ifdef MFEM_THREAD_SAFE
1164  DenseMatrix dshape(test_nd, dim);
1165  DenseMatrix dshapedxt(test_nd, dim);
1166  DenseMatrix vshape(trial_nd, dim);
1167  DenseMatrix invdfdx(dim);
1168 #else
1169  dshape.SetSize(test_nd, dim);
1170  dshapedxt.SetSize(test_nd, dim);
1171  vshape.SetSize(trial_nd, dim);
1172  invdfdx.SetSize(dim);
1173 #endif
1174 
1175  elmat.SetSize(test_nd, trial_nd);
1176 
1177  const IntegrationRule *ir = IntRule;
1178  if (ir == NULL)
1179  {
1180  // The integrand on the reference element is:
1181  // -( Q/det(J) ) u_hat^T adj(J) adj(J)^T grad_hat(v_hat).
1182  //
1183  // For Trans in (P_k)^d, v_hat in P_l, u_hat in ND_m, and dim=sdim=d>=1
1184  // - J_{ij} is in P_{k-1}, so adj(J)_{ij} is in P_{(d-1)*(k-1)}
1185  // - so adj(J)^T grad_hat(v_hat) is in (P_{(d-1)*(k-1)+(l-1)})^d
1186  // - u_hat is in (P_m)^d
1187  // - adj(J)^T u_hat is in (P_{(d-1)*(k-1)+m})^d
1188  // - and u_hat^T adj(J) adj(J)^T grad_hat(v_hat) is in P_n with
1189  // n = 2*(d-1)*(k-1)+(l-1)+m
1190  //
1191  // For Trans in (Q_k)^d, v_hat in Q_l, u_hat in ND_m, and dim=sdim=d>1
1192  // - J_{i*}, J's i-th row, is in ( Q_{k-1,k,k}, Q_{k,k-1,k}, Q_{k,k,k-1} )
1193  // - adj(J)_{*j} is in ( Q_{s,s-1,s-1}, Q_{s-1,s,s-1}, Q_{s-1,s-1,s} )
1194  // with s = (d-1)*k
1195  // - adj(J)^T grad_hat(v_hat) is in Q_{(d-1)*k+(l-1)}
1196  // - u_hat is in ( Q_{m-1,m,m}, Q_{m,m-1,m}, Q_{m,m,m-1} )
1197  // - adj(J)^T u_hat is in Q_{(d-1)*k+(m-1)}
1198  // - and u_hat^T adj(J) adj(J)^T grad_hat(v_hat) is in Q_n with
1199  // n = 2*(d-1)*k+(l-1)+(m-1)
1200  //
1201  // In the next formula we use the expressions for n with k=1, which means
1202  // that the term Q/det(J) is disregard:
1203  int ir_order = (trial_fe.Space() == FunctionSpace::Pk) ?
1204  (trial_fe.GetOrder() + test_fe.GetOrder() - 1) :
1205  (trial_fe.GetOrder() + test_fe.GetOrder() + 2*(dim-2));
1206  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
1207  }
1208 
1209  elmat = 0.0;
1210  for (i = 0; i < ir->GetNPoints(); i++)
1211  {
1212  const IntegrationPoint &ip = ir->IntPoint(i);
1213  test_fe.CalcDShape(ip, dshape);
1214 
1215  Trans.SetIntPoint(&ip);
1216  CalcAdjugate(Trans.Jacobian(), invdfdx);
1217  Mult(dshape, invdfdx, dshapedxt);
1218 
1219  trial_fe.CalcVShape(Trans, vshape);
1220 
1221  double w = ip.weight;
1222 
1223  if (Q)
1224  {
1225  w *= Q->Eval(Trans, ip);
1226  }
1227  dshapedxt *= -w;
1228 
1229  AddMultABt(dshapedxt, vshape, elmat);
1230  }
1231 }
1232 
1233 void VectorFECurlIntegrator::AssembleElementMatrix2(
1234  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1235  ElementTransformation &Trans, DenseMatrix &elmat)
1236 {
1237  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1238  int dim = trial_fe.GetDim();
1239  int dimc = (dim == 3) ? 3 : 1;
1240 
1241  MFEM_ASSERT(trial_fe.GetMapType() == mfem::FiniteElement::H_CURL ||
1243  "At least one of the finite elements must be in H(Curl)");
1244 
1245  int curl_nd, vec_nd;
1246  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1247  {
1248  curl_nd = trial_nd;
1249  vec_nd = test_nd;
1250  }
1251  else
1252  {
1253  curl_nd = test_nd;
1254  vec_nd = trial_nd;
1255  }
1256 
1257 #ifdef MFEM_THREAD_SAFE
1258  DenseMatrix curlshapeTrial(curl_nd, dimc);
1259  DenseMatrix curlshapeTrial_dFT(curl_nd, dimc);
1260  DenseMatrix vshapeTest(vec_nd, dimc);
1261 #else
1262  curlshapeTrial.SetSize(curl_nd, dimc);
1263  curlshapeTrial_dFT.SetSize(curl_nd, dimc);
1264  vshapeTest.SetSize(vec_nd, dimc);
1265 #endif
1266  Vector shapeTest(vshapeTest.GetData(), vec_nd);
1267 
1268  elmat.SetSize(test_nd, trial_nd);
1269 
1270  const IntegrationRule *ir = IntRule;
1271  if (ir == NULL)
1272  {
1273  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
1274  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1275  }
1276 
1277  elmat = 0.0;
1278  for (i = 0; i < ir->GetNPoints(); i++)
1279  {
1280  const IntegrationPoint &ip = ir->IntPoint(i);
1281 
1282  Trans.SetIntPoint(&ip);
1283  if (dim == 3)
1284  {
1285  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1286  {
1287  trial_fe.CalcCurlShape(ip, curlshapeTrial);
1288  test_fe.CalcVShape(Trans, vshapeTest);
1289  }
1290  else
1291  {
1292  test_fe.CalcCurlShape(ip, curlshapeTrial);
1293  trial_fe.CalcVShape(Trans, vshapeTest);
1294  }
1295  MultABt(curlshapeTrial, Trans.Jacobian(), curlshapeTrial_dFT);
1296  }
1297  else
1298  {
1299  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1300  {
1301  trial_fe.CalcCurlShape(ip, curlshapeTrial_dFT);
1302  test_fe.CalcShape(ip, shapeTest);
1303  }
1304  else
1305  {
1306  test_fe.CalcCurlShape(ip, curlshapeTrial_dFT);
1307  trial_fe.CalcShape(ip, shapeTest);
1308  }
1309  }
1310 
1311  double w = ip.weight;
1312 
1313  if (Q)
1314  {
1315  w *= Q->Eval(Trans, ip);
1316  }
1317  // Note: shapeTest points to the same data as vshapeTest
1318  vshapeTest *= w;
1319  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1320  {
1321  AddMultABt(vshapeTest, curlshapeTrial_dFT, elmat);
1322  }
1323  else
1324  {
1325  AddMultABt(curlshapeTrial_dFT, vshapeTest, elmat);
1326  }
1327  }
1328 }
1329 
1330 void DerivativeIntegrator::AssembleElementMatrix2 (
1331  const FiniteElement &trial_fe,
1332  const FiniteElement &test_fe,
1333  ElementTransformation &Trans,
1334  DenseMatrix &elmat)
1335 {
1336  int dim = trial_fe.GetDim();
1337  int trial_nd = trial_fe.GetDof();
1338  int test_nd = test_fe.GetDof();
1339 
1340  int i, l;
1341  double det;
1342 
1343  elmat.SetSize (test_nd,trial_nd);
1344  dshape.SetSize (trial_nd,dim);
1345  dshapedxt.SetSize(trial_nd,dim);
1346  dshapedxi.SetSize(trial_nd);
1347  invdfdx.SetSize(dim);
1348  shape.SetSize (test_nd);
1349 
1350  const IntegrationRule *ir = IntRule;
1351  if (ir == NULL)
1352  {
1353  int order;
1354  if (trial_fe.Space() == FunctionSpace::Pk)
1355  {
1356  order = trial_fe.GetOrder() + test_fe.GetOrder() - 1;
1357  }
1358  else
1359  {
1360  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim;
1361  }
1362 
1363  if (trial_fe.Space() == FunctionSpace::rQk)
1364  {
1365  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1366  }
1367  else
1368  {
1369  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1370  }
1371  }
1372 
1373  elmat = 0.0;
1374  for (i = 0; i < ir->GetNPoints(); i++)
1375  {
1376  const IntegrationPoint &ip = ir->IntPoint(i);
1377 
1378  trial_fe.CalcDShape(ip, dshape);
1379 
1380  Trans.SetIntPoint (&ip);
1381  CalcInverse (Trans.Jacobian(), invdfdx);
1382  det = Trans.Weight();
1383  Mult (dshape, invdfdx, dshapedxt);
1384 
1385  test_fe.CalcShape(ip, shape);
1386 
1387  for (l = 0; l < trial_nd; l++)
1388  {
1389  dshapedxi(l) = dshapedxt(l,xi);
1390  }
1391 
1392  shape *= Q.Eval(Trans,ip) * det * ip.weight;
1393  AddMultVWt (shape, dshapedxi, elmat);
1394  }
1395 }
1396 
1397 void CurlCurlIntegrator::AssembleElementMatrix
1399  DenseMatrix &elmat )
1400 {
1401  int nd = el.GetDof();
1402  int dim = el.GetDim();
1403  int dimc = (dim == 3) ? 3 : 1;
1404  double w;
1405 
1406 #ifdef MFEM_THREAD_SAFE
1407  DenseMatrix curlshape(nd,dimc), curlshape_dFt(nd,dimc);
1408 #else
1409  curlshape.SetSize(nd,dimc);
1410  curlshape_dFt.SetSize(nd,dimc);
1411 #endif
1412  elmat.SetSize(nd);
1413 
1414  const IntegrationRule *ir = IntRule;
1415  if (ir == NULL)
1416  {
1417  int order;
1418  if (el.Space() == FunctionSpace::Pk)
1419  {
1420  order = 2*el.GetOrder() - 2;
1421  }
1422  else
1423  {
1424  order = 2*el.GetOrder();
1425  }
1426 
1427  ir = &IntRules.Get(el.GetGeomType(), order);
1428  }
1429 
1430  elmat = 0.0;
1431  for (int i = 0; i < ir->GetNPoints(); i++)
1432  {
1433  const IntegrationPoint &ip = ir->IntPoint(i);
1434 
1435  Trans.SetIntPoint (&ip);
1436 
1437  w = ip.weight / Trans.Weight();
1438 
1439  if ( dim == 3 )
1440  {
1441  el.CalcCurlShape(ip, curlshape);
1442  MultABt(curlshape, Trans.Jacobian(), curlshape_dFt);
1443  }
1444  else
1445  {
1446  el.CalcCurlShape(ip, curlshape_dFt);
1447  }
1448 
1449  if (Q)
1450  {
1451  w *= Q->Eval(Trans, ip);
1452  }
1453 
1454  AddMult_a_AAt(w, curlshape_dFt, elmat);
1455  }
1456 }
1457 
1458 void CurlCurlIntegrator
1459 ::ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans,
1460  Vector &u, const FiniteElement &fluxelem, Vector &flux,
1461  int with_coef)
1462 {
1463 #ifdef MFEM_THREAD_SAFE
1464  DenseMatrix projcurl;
1465 #endif
1466 
1467  fluxelem.ProjectCurl(el, Trans, projcurl);
1468 
1469  flux.SetSize(projcurl.Height());
1470  projcurl.Mult(u, flux);
1471 
1472  // TODO: Q, wcoef?
1473 }
1474 
1475 double CurlCurlIntegrator::ComputeFluxEnergy(const FiniteElement &fluxelem,
1476  ElementTransformation &Trans,
1477  Vector &flux, Vector *d_energy)
1478 {
1479  int nd = fluxelem.GetDof();
1480  int dim = fluxelem.GetDim();
1481 
1482 #ifdef MFEM_THREAD_SAFE
1483  DenseMatrix vshape;
1484 #endif
1485  vshape.SetSize(nd, dim);
1486  pointflux.SetSize(dim);
1487  if (d_energy) { vec.SetSize(dim); }
1488 
1489  int order = 2 * fluxelem.GetOrder(); // <--
1490  const IntegrationRule &ir = IntRules.Get(fluxelem.GetGeomType(), order);
1491 
1492  double energy = 0.0;
1493  if (d_energy) { *d_energy = 0.0; }
1494 
1495  Vector* pfluxes = NULL;
1496  if (d_energy)
1497  {
1498  pfluxes = new Vector[ir.GetNPoints()];
1499  }
1500 
1501  for (int i = 0; i < ir.GetNPoints(); i++)
1502  {
1503  const IntegrationPoint &ip = ir.IntPoint(i);
1504  Trans.SetIntPoint(&ip);
1505 
1506  fluxelem.CalcVShape(Trans, vshape);
1507  // fluxelem.CalcVShape(ip, vshape);
1508  vshape.MultTranspose(flux, pointflux);
1509 
1510  double w = Trans.Weight() * ip.weight;
1511 
1512  double e = w * (pointflux * pointflux);
1513 
1514  if (Q)
1515  {
1516  // TODO
1517  }
1518 
1519  energy += e;
1520 
1521 #if ANISO_EXPERIMENTAL
1522  if (d_energy)
1523  {
1524  pfluxes[i].SetSize(dim);
1525  Trans.Jacobian().MultTranspose(pointflux, pfluxes[i]);
1526 
1527  /*
1528  DenseMatrix Jadj(dim, dim);
1529  CalcAdjugate(Trans.Jacobian(), Jadj);
1530  pfluxes[i].SetSize(dim);
1531  Jadj.Mult(pointflux, pfluxes[i]);
1532  */
1533 
1534  // pfluxes[i] = pointflux;
1535  }
1536 #endif
1537  }
1538 
1539  if (d_energy)
1540  {
1541 #if ANISO_EXPERIMENTAL
1542  *d_energy = 0.0;
1543  Vector tmp;
1544 
1545  int n = (int) round(pow(ir.GetNPoints(), 1.0/3.0));
1546  MFEM_ASSERT(n*n*n == ir.GetNPoints(), "");
1547 
1548  // hack: get total variation of 'pointflux' in the x,y,z directions
1549  for (int k = 0; k < n; k++)
1550  for (int l = 0; l < n; l++)
1551  for (int m = 0; m < n; m++)
1552  {
1553  Vector &vec = pfluxes[(k*n + l)*n + m];
1554  if (m > 0)
1555  {
1556  tmp = vec; tmp -= pfluxes[(k*n + l)*n + (m-1)];
1557  (*d_energy)[0] += (tmp * tmp);
1558  }
1559  if (l > 0)
1560  {
1561  tmp = vec; tmp -= pfluxes[(k*n + (l-1))*n + m];
1562  (*d_energy)[1] += (tmp * tmp);
1563  }
1564  if (k > 0)
1565  {
1566  tmp = vec; tmp -= pfluxes[((k-1)*n + l)*n + m];
1567  (*d_energy)[2] += (tmp * tmp);
1568  }
1569  }
1570 #else
1571  *d_energy = 1.0;
1572 #endif
1573 
1574  delete [] pfluxes;
1575  }
1576 
1577  return energy;
1578 }
1579 
1580 void VectorCurlCurlIntegrator::AssembleElementMatrix(
1581  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1582 {
1583  int dim = el.GetDim();
1584  int dof = el.GetDof();
1585  int cld = (dim*(dim-1))/2;
1586 
1587 #ifdef MFEM_THREAD_SAFE
1588  DenseMatrix dshape_hat(dof, dim), dshape(dof, dim);
1589  DenseMatrix curlshape(dim*dof, cld), Jadj(dim);
1590 #else
1591  dshape_hat.SetSize(dof, dim);
1592  dshape.SetSize(dof, dim);
1593  curlshape.SetSize(dim*dof, cld);
1594  Jadj.SetSize(dim);
1595 #endif
1596 
1597  const IntegrationRule *ir = IntRule;
1598  if (ir == NULL)
1599  {
1600  // use the same integration rule as diffusion
1601  int order = 2 * Trans.OrderGrad(&el);
1602  ir = &IntRules.Get(el.GetGeomType(), order);
1603  }
1604 
1605  elmat = 0.0;
1606  for (int i = 0; i < ir->GetNPoints(); i++)
1607  {
1608  const IntegrationPoint &ip = ir->IntPoint(i);
1609  el.CalcDShape(ip, dshape_hat);
1610 
1611  Trans.SetIntPoint(&ip);
1612  CalcAdjugate(Trans.Jacobian(), Jadj);
1613  double w = ip.weight / Trans.Weight();
1614 
1615  Mult(dshape_hat, Jadj, dshape);
1616  dshape.GradToCurl(curlshape);
1617 
1618  if (Q)
1619  {
1620  w *= Q->Eval(Trans, ip);
1621  }
1622 
1623  AddMult_a_AAt(w, curlshape, elmat);
1624  }
1625 }
1626 
1627 double VectorCurlCurlIntegrator::GetElementEnergy(
1628  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
1629 {
1630  int dim = el.GetDim();
1631  int dof = el.GetDof();
1632 
1633 #ifdef MFEM_THREAD_SAFE
1634  DenseMatrix dshape_hat(dof, dim), Jadj(dim), grad_hat(dim), grad(dim);
1635 #else
1636  dshape_hat.SetSize(dof, dim);
1637 
1638  Jadj.SetSize(dim);
1639  grad_hat.SetSize(dim);
1640  grad.SetSize(dim);
1641 #endif
1642  DenseMatrix elfun_mat(elfun.GetData(), dof, dim);
1643 
1644  const IntegrationRule *ir = IntRule;
1645  if (ir == NULL)
1646  {
1647  // use the same integration rule as diffusion
1648  int order = 2 * Tr.OrderGrad(&el);
1649  ir = &IntRules.Get(el.GetGeomType(), order);
1650  }
1651 
1652  double energy = 0.;
1653  for (int i = 0; i < ir->GetNPoints(); i++)
1654  {
1655  const IntegrationPoint &ip = ir->IntPoint(i);
1656  el.CalcDShape(ip, dshape_hat);
1657 
1658  MultAtB(elfun_mat, dshape_hat, grad_hat);
1659 
1660  Tr.SetIntPoint(&ip);
1661  CalcAdjugate(Tr.Jacobian(), Jadj);
1662  double w = ip.weight / Tr.Weight();
1663 
1664  Mult(grad_hat, Jadj, grad);
1665 
1666  if (dim == 2)
1667  {
1668  double curl = grad(0,1) - grad(1,0);
1669  w *= curl * curl;
1670  }
1671  else
1672  {
1673  double curl_x = grad(2,1) - grad(1,2);
1674  double curl_y = grad(0,2) - grad(2,0);
1675  double curl_z = grad(1,0) - grad(0,1);
1676  w *= curl_x * curl_x + curl_y * curl_y + curl_z * curl_z;
1677  }
1678 
1679  if (Q)
1680  {
1681  w *= Q->Eval(Tr, ip);
1682  }
1683 
1684  energy += w;
1685  }
1686 
1687  elfun_mat.ClearExternalData();
1688 
1689  return 0.5 * energy;
1690 }
1691 
1692 
1693 void VectorFEMassIntegrator::AssembleElementMatrix(
1694  const FiniteElement &el,
1695  ElementTransformation &Trans,
1696  DenseMatrix &elmat)
1697 {
1698  int dof = el.GetDof();
1699  int spaceDim = Trans.GetSpaceDim();
1700 
1701  double w;
1702 
1703 #ifdef MFEM_THREAD_SAFE
1704  Vector D(VQ ? VQ->GetVDim() : 0);
1705  DenseMatrix trial_vshape(dof, spaceDim);
1706  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1707 #else
1708  trial_vshape.SetSize(dof, spaceDim);
1709  D.SetSize(VQ ? VQ->GetVDim() : 0);
1710  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1711 #endif
1712  DenseMatrix tmp(trial_vshape.Height(), K.Width());
1713 
1714  elmat.SetSize(dof);
1715  elmat = 0.0;
1716 
1717  const IntegrationRule *ir = IntRule;
1718  if (ir == NULL)
1719  {
1720  // int order = 2 * el.GetOrder();
1721  int order = Trans.OrderW() + 2 * el.GetOrder();
1722  ir = &IntRules.Get(el.GetGeomType(), order);
1723  }
1724 
1725  for (int i = 0; i < ir->GetNPoints(); i++)
1726  {
1727  const IntegrationPoint &ip = ir->IntPoint(i);
1728 
1729  Trans.SetIntPoint (&ip);
1730 
1731  el.CalcVShape(Trans, trial_vshape);
1732 
1733  w = ip.weight * Trans.Weight();
1734  if (MQ)
1735  {
1736  MQ->Eval(K, Trans, ip);
1737  K *= w;
1738  Mult(trial_vshape,K,tmp);
1739  AddMultABt(tmp,trial_vshape,elmat);
1740  }
1741  else if (VQ)
1742  {
1743  VQ->Eval(D, Trans, ip);
1744  D *= w;
1745  AddMultADAt(trial_vshape, D, elmat);
1746  }
1747  else
1748  {
1749  if (Q)
1750  {
1751  w *= Q -> Eval (Trans, ip);
1752  }
1753  AddMult_a_AAt (w, trial_vshape, elmat);
1754  }
1755  }
1756 }
1757 
1758 void VectorFEMassIntegrator::AssembleElementMatrix2(
1759  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1760  ElementTransformation &Trans, DenseMatrix &elmat)
1761 {
1762  if ( test_fe.GetRangeType() == FiniteElement::SCALAR && VQ )
1763  {
1764  // assume test_fe is scalar FE and trial_fe is vector FE
1765  int dim = test_fe.GetDim();
1766  int trial_dof = trial_fe.GetDof();
1767  int test_dof = test_fe.GetDof();
1768  double w;
1769 
1770  if (MQ)
1771  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1772  " is not implemented for tensor materials");
1773 
1774 #ifdef MFEM_THREAD_SAFE
1775  DenseMatrix trial_vshape(trial_dof, dim);
1776  Vector shape(test_dof);
1777  Vector D(dim);
1778 #else
1779  trial_vshape.SetSize(trial_dof, dim);
1780  shape.SetSize(test_dof);
1781  D.SetSize(dim);
1782 #endif
1783 
1784  elmat.SetSize (test_dof, trial_dof);
1785 
1786  const IntegrationRule *ir = IntRule;
1787  if (ir == NULL)
1788  {
1789  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1790  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1791  }
1792 
1793  elmat = 0.0;
1794  for (int i = 0; i < ir->GetNPoints(); i++)
1795  {
1796  const IntegrationPoint &ip = ir->IntPoint(i);
1797 
1798  Trans.SetIntPoint (&ip);
1799 
1800  trial_fe.CalcVShape(Trans, trial_vshape);
1801  test_fe.CalcShape(ip, shape);
1802 
1803  w = ip.weight * Trans.Weight();
1804  VQ->Eval(D, Trans, ip);
1805  D *= w;
1806 
1807  for (int d = 0; d < dim; d++)
1808  {
1809  for (int j = 0; j < test_dof; j++)
1810  {
1811  for (int k = 0; k < trial_dof; k++)
1812  {
1813  elmat(j, k) += D[d] * shape(j) * trial_vshape(k, d);
1814  }
1815  }
1816  }
1817  }
1818  }
1819  else if ( test_fe.GetRangeType() == FiniteElement::SCALAR )
1820  {
1821  // assume test_fe is scalar FE and trial_fe is vector FE
1822  int dim = test_fe.GetDim();
1823  int trial_dof = trial_fe.GetDof();
1824  int test_dof = test_fe.GetDof();
1825  double w;
1826 
1827  if (VQ || MQ)
1828  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1829  " is not implemented for vector/tensor permeability");
1830 
1831 #ifdef MFEM_THREAD_SAFE
1832  DenseMatrix trial_vshape(trial_dof, dim);
1833  Vector shape(test_dof);
1834 #else
1835  trial_vshape.SetSize(trial_dof, dim);
1836  shape.SetSize(test_dof);
1837 #endif
1838 
1839  elmat.SetSize (dim*test_dof, trial_dof);
1840 
1841  const IntegrationRule *ir = IntRule;
1842  if (ir == NULL)
1843  {
1844  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1845  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1846  }
1847 
1848  elmat = 0.0;
1849  for (int i = 0; i < ir->GetNPoints(); i++)
1850  {
1851  const IntegrationPoint &ip = ir->IntPoint(i);
1852 
1853  Trans.SetIntPoint (&ip);
1854 
1855  trial_fe.CalcVShape(Trans, trial_vshape);
1856  test_fe.CalcShape(ip, shape);
1857 
1858  w = ip.weight * Trans.Weight();
1859  if (Q)
1860  {
1861  w *= Q -> Eval (Trans, ip);
1862  }
1863 
1864  for (int d = 0; d < dim; d++)
1865  {
1866  for (int j = 0; j < test_dof; j++)
1867  {
1868  for (int k = 0; k < trial_dof; k++)
1869  {
1870  elmat(d * test_dof + j, k) += w * shape(j) * trial_vshape(k, d);
1871  }
1872  }
1873  }
1874  }
1875  }
1876  else
1877  {
1878  // assume both test_fe and trial_fe are vector FE
1879  int dim = test_fe.GetDim();
1880  int trial_dof = trial_fe.GetDof();
1881  int test_dof = test_fe.GetDof();
1882  double w;
1883 
1884  if (VQ || MQ)
1885  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1886  " is not implemented for vector/tensor permeability");
1887 
1888 #ifdef MFEM_THREAD_SAFE
1889  DenseMatrix trial_vshape(trial_dof, dim);
1890  DenseMatrix test_vshape(test_dof,dim);
1891 #else
1892  trial_vshape.SetSize(trial_dof, dim);
1893  test_vshape.SetSize(test_dof,dim);
1894 #endif
1895 
1896  elmat.SetSize (test_dof, trial_dof);
1897 
1898  const IntegrationRule *ir = IntRule;
1899  if (ir == NULL)
1900  {
1901  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1902  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1903  }
1904 
1905  elmat = 0.0;
1906  for (int i = 0; i < ir->GetNPoints(); i++)
1907  {
1908  const IntegrationPoint &ip = ir->IntPoint(i);
1909 
1910  Trans.SetIntPoint (&ip);
1911 
1912  trial_fe.CalcVShape(Trans, trial_vshape);
1913  test_fe.CalcVShape(Trans, test_vshape);
1914 
1915  w = ip.weight * Trans.Weight();
1916  if (Q)
1917  {
1918  w *= Q -> Eval (Trans, ip);
1919  }
1920 
1921  for (int d = 0; d < dim; d++)
1922  {
1923  for (int j = 0; j < test_dof; j++)
1924  {
1925  for (int k = 0; k < trial_dof; k++)
1926  {
1927  elmat(j, k) += w * test_vshape(j, d) * trial_vshape(k, d);
1928  }
1929  }
1930  }
1931  }
1932  }
1933 }
1934 
1935 void VectorDivergenceIntegrator::AssembleElementMatrix2(
1936  const FiniteElement &trial_fe,
1937  const FiniteElement &test_fe,
1938  ElementTransformation &Trans,
1939  DenseMatrix &elmat)
1940 {
1941  int dim = trial_fe.GetDim();
1942  int trial_dof = trial_fe.GetDof();
1943  int test_dof = test_fe.GetDof();
1944  double c;
1945 
1946  dshape.SetSize (trial_dof, dim);
1947  gshape.SetSize (trial_dof, dim);
1948  Jadj.SetSize (dim);
1949  divshape.SetSize (dim*trial_dof);
1950  shape.SetSize (test_dof);
1951 
1952  elmat.SetSize (test_dof, dim*trial_dof);
1953 
1954  const IntegrationRule *ir = IntRule;
1955  if (ir == NULL)
1956  {
1957  int order = Trans.OrderGrad(&trial_fe) + test_fe.GetOrder();
1958  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1959  }
1960 
1961  elmat = 0.0;
1962 
1963  for (int i = 0; i < ir -> GetNPoints(); i++)
1964  {
1965  const IntegrationPoint &ip = ir->IntPoint(i);
1966 
1967  trial_fe.CalcDShape (ip, dshape);
1968  test_fe.CalcShape (ip, shape);
1969 
1970  Trans.SetIntPoint (&ip);
1971  CalcAdjugate(Trans.Jacobian(), Jadj);
1972 
1973  Mult (dshape, Jadj, gshape);
1974 
1975  gshape.GradToDiv (divshape);
1976 
1977  c = ip.weight;
1978  if (Q)
1979  {
1980  c *= Q -> Eval (Trans, ip);
1981  }
1982 
1983  // elmat += c * shape * divshape ^ t
1984  shape *= c;
1985  AddMultVWt (shape, divshape, elmat);
1986  }
1987 }
1988 
1989 
1990 void DivDivIntegrator::AssembleElementMatrix(
1991  const FiniteElement &el,
1992  ElementTransformation &Trans,
1993  DenseMatrix &elmat)
1994 {
1995  int dof = el.GetDof();
1996  double c;
1997 
1998 #ifdef MFEM_THREAD_SAFE
1999  Vector divshape(dof);
2000 #else
2001  divshape.SetSize(dof);
2002 #endif
2003  elmat.SetSize(dof);
2004 
2005  const IntegrationRule *ir = IntRule;
2006  if (ir == NULL)
2007  {
2008  int order = 2 * el.GetOrder() - 2; // <--- OK for RTk
2009  ir = &IntRules.Get(el.GetGeomType(), order);
2010  }
2011 
2012  elmat = 0.0;
2013 
2014  for (int i = 0; i < ir -> GetNPoints(); i++)
2015  {
2016  const IntegrationPoint &ip = ir->IntPoint(i);
2017 
2018  el.CalcDivShape (ip, divshape);
2019 
2020  Trans.SetIntPoint (&ip);
2021  c = ip.weight / Trans.Weight();
2022 
2023  if (Q)
2024  {
2025  c *= Q -> Eval (Trans, ip);
2026  }
2027 
2028  // elmat += c * divshape * divshape ^ t
2029  AddMult_a_VVt (c, divshape, elmat);
2030  }
2031 }
2032 
2033 
2034 void VectorDiffusionIntegrator::AssembleElementMatrix(
2035  const FiniteElement &el,
2036  ElementTransformation &Trans,
2037  DenseMatrix &elmat)
2038 {
2039  int dim = el.GetDim();
2040  int dof = el.GetDof();
2041 
2042  double norm;
2043 
2044  elmat.SetSize (dim * dof);
2045 
2046  Jinv. SetSize (dim);
2047  dshape.SetSize (dof, dim);
2048  gshape.SetSize (dof, dim);
2049  pelmat.SetSize (dof);
2050 
2051  const IntegrationRule *ir = IntRule;
2052  if (ir == NULL)
2053  {
2054  // integrand is rational function if det(J) is not constant
2055  int order = 2 * Trans.OrderGrad(&el); // order of the numerator
2056  if (el.Space() == FunctionSpace::rQk)
2057  {
2058  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
2059  }
2060  else
2061  {
2062  ir = &IntRules.Get(el.GetGeomType(), order);
2063  }
2064  }
2065 
2066  elmat = 0.0;
2067 
2068  for (int i = 0; i < ir -> GetNPoints(); i++)
2069  {
2070  const IntegrationPoint &ip = ir->IntPoint(i);
2071 
2072  el.CalcDShape (ip, dshape);
2073 
2074  Trans.SetIntPoint (&ip);
2075  norm = ip.weight * Trans.Weight();
2076  CalcInverse (Trans.Jacobian(), Jinv);
2077 
2078  Mult (dshape, Jinv, gshape);
2079 
2080  MultAAt (gshape, pelmat);
2081 
2082  if (Q)
2083  {
2084  norm *= Q -> Eval (Trans, ip);
2085  }
2086 
2087  pelmat *= norm;
2088 
2089  for (int d = 0; d < dim; d++)
2090  {
2091  for (int k = 0; k < dof; k++)
2092  for (int l = 0; l < dof; l++)
2093  {
2094  elmat (dof*d+k, dof*d+l) += pelmat (k, l);
2095  }
2096  }
2097  }
2098 }
2099 
2100 void VectorDiffusionIntegrator::AssembleElementVector(
2101  const FiniteElement &el, ElementTransformation &Tr,
2102  const Vector &elfun, Vector &elvect)
2103 {
2104  int dim = el.GetDim(); // assuming vector_dim == reference_dim
2105  int dof = el.GetDof();
2106  double w;
2107 
2108  Jinv.SetSize(dim);
2109  dshape.SetSize(dof, dim);
2110  pelmat.SetSize(dim);
2111  gshape.SetSize(dim);
2112 
2113  elvect.SetSize(dim*dof);
2114  DenseMatrix mat_in(elfun.GetData(), dof, dim);
2115  DenseMatrix mat_out(elvect.GetData(), dof, dim);
2116 
2117  const IntegrationRule *ir = IntRule;
2118  if (ir == NULL)
2119  {
2120  // integrant is rational function if det(J) is not constant
2121  int order = 2 * Tr.OrderGrad(&el); // order of the numerator
2122  ir = (el.Space() == FunctionSpace::rQk) ?
2123  &RefinedIntRules.Get(el.GetGeomType(), order) :
2124  &IntRules.Get(el.GetGeomType(), order);
2125  }
2126 
2127  elvect = 0.0;
2128  for (int i = 0; i < ir->GetNPoints(); i++)
2129  {
2130  const IntegrationPoint &ip = ir->IntPoint(i);
2131 
2132  Tr.SetIntPoint(&ip);
2133  CalcAdjugate(Tr.Jacobian(), Jinv);
2134  w = ip.weight / Tr.Weight();
2135  if (Q)
2136  {
2137  w *= Q->Eval(Tr, ip);
2138  }
2139  MultAAt(Jinv, gshape);
2140  gshape *= w;
2141 
2142  el.CalcDShape(ip, dshape);
2143 
2144  MultAtB(mat_in, dshape, pelmat);
2145  MultABt(pelmat, gshape, Jinv);
2146  AddMultABt(dshape, Jinv, mat_out);
2147  }
2148 }
2149 
2150 
2151 void ElasticityIntegrator::AssembleElementMatrix(
2152  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
2153 {
2154  int dof = el.GetDof();
2155  int dim = el.GetDim();
2156  double w, L, M;
2157 
2158 #ifdef MFEM_THREAD_SAFE
2159  DenseMatrix dshape(dof, dim), Jinv(dim), gshape(dof, dim), pelmat(dof);
2160  Vector divshape(dim*dof);
2161 #else
2162  Jinv.SetSize(dim);
2163  dshape.SetSize(dof, dim);
2164  gshape.SetSize(dof, dim);
2165  pelmat.SetSize(dof);
2166  divshape.SetSize(dim*dof);
2167 #endif
2168 
2169  elmat.SetSize(dof * dim);
2170 
2171  const IntegrationRule *ir = IntRule;
2172  if (ir == NULL)
2173  {
2174  int order = 2 * Trans.OrderGrad(&el); // correct order?
2175  ir = &IntRules.Get(el.GetGeomType(), order);
2176  }
2177 
2178  elmat = 0.0;
2179 
2180  for (int i = 0; i < ir -> GetNPoints(); i++)
2181  {
2182  const IntegrationPoint &ip = ir->IntPoint(i);
2183 
2184  el.CalcDShape(ip, dshape);
2185 
2186  Trans.SetIntPoint(&ip);
2187  w = ip.weight * Trans.Weight();
2188  CalcInverse(Trans.Jacobian(), Jinv);
2189  Mult(dshape, Jinv, gshape);
2190  MultAAt(gshape, pelmat);
2191  gshape.GradToDiv (divshape);
2192 
2193  M = mu->Eval(Trans, ip);
2194  if (lambda)
2195  {
2196  L = lambda->Eval(Trans, ip);
2197  }
2198  else
2199  {
2200  L = q_lambda * M;
2201  M = q_mu * M;
2202  }
2203 
2204  if (L != 0.0)
2205  {
2206  AddMult_a_VVt(L * w, divshape, elmat);
2207  }
2208 
2209  if (M != 0.0)
2210  {
2211  for (int d = 0; d < dim; d++)
2212  {
2213  for (int k = 0; k < dof; k++)
2214  for (int l = 0; l < dof; l++)
2215  {
2216  elmat (dof*d+k, dof*d+l) += (M * w) * pelmat(k, l);
2217  }
2218  }
2219  for (int i = 0; i < dim; i++)
2220  for (int j = 0; j < dim; j++)
2221  {
2222  for (int k = 0; k < dof; k++)
2223  for (int l = 0; l < dof; l++)
2224  elmat(dof*i+k, dof*j+l) +=
2225  (M * w) * gshape(k, j) * gshape(l, i);
2226  // + (L * w) * gshape(k, i) * gshape(l, j)
2227  }
2228  }
2229  }
2230 }
2231 
2232 void DGTraceIntegrator::AssembleFaceMatrix(const FiniteElement &el1,
2233  const FiniteElement &el2,
2235  DenseMatrix &elmat)
2236 {
2237  int dim, ndof1, ndof2;
2238 
2239  double un, a, b, w;
2240 
2241  dim = el1.GetDim();
2242  ndof1 = el1.GetDof();
2243  Vector vu(dim), nor(dim);
2244 
2245  if (Trans.Elem2No >= 0)
2246  {
2247  ndof2 = el2.GetDof();
2248  }
2249  else
2250  {
2251  ndof2 = 0;
2252  }
2253 
2254  shape1.SetSize(ndof1);
2255  shape2.SetSize(ndof2);
2256  elmat.SetSize(ndof1 + ndof2);
2257  elmat = 0.0;
2258 
2259  const IntegrationRule *ir = IntRule;
2260  if (ir == NULL)
2261  {
2262  int order;
2263  // Assuming order(u)==order(mesh)
2264  if (Trans.Elem2No >= 0)
2265  order = (min(Trans.Elem1->OrderW(), Trans.Elem2->OrderW()) +
2266  2*max(el1.GetOrder(), el2.GetOrder()));
2267  else
2268  {
2269  order = Trans.Elem1->OrderW() + 2*el1.GetOrder();
2270  }
2271  if (el1.Space() == FunctionSpace::Pk)
2272  {
2273  order++;
2274  }
2275  ir = &IntRules.Get(Trans.FaceGeom, order);
2276  }
2277 
2278  for (int p = 0; p < ir->GetNPoints(); p++)
2279  {
2280  const IntegrationPoint &ip = ir->IntPoint(p);
2281  IntegrationPoint eip1, eip2;
2282  Trans.Loc1.Transform(ip, eip1);
2283  if (ndof2)
2284  {
2285  Trans.Loc2.Transform(ip, eip2);
2286  }
2287  el1.CalcShape(eip1, shape1);
2288 
2289  Trans.Face->SetIntPoint(&ip);
2290  Trans.Elem1->SetIntPoint(&eip1);
2291 
2292  u->Eval(vu, *Trans.Elem1, eip1);
2293 
2294  if (dim == 1)
2295  {
2296  nor(0) = 2*eip1.x - 1.0;
2297  }
2298  else
2299  {
2300  CalcOrtho(Trans.Face->Jacobian(), nor);
2301  }
2302 
2303  un = vu * nor;
2304  a = 0.5 * alpha * un;
2305  b = beta * fabs(un);
2306  // note: if |alpha/2|==|beta| then |a|==|b|, i.e. (a==b) or (a==-b)
2307  // and therefore two blocks in the element matrix contribution
2308  // (from the current quadrature point) are 0
2309 
2310  if (rho)
2311  {
2312  double rho_p;
2313  if (un >= 0.0 && ndof2)
2314  {
2315  Trans.Elem2->SetIntPoint(&eip2);
2316  rho_p = rho->Eval(*Trans.Elem2, eip2);
2317  }
2318  else
2319  {
2320  rho_p = rho->Eval(*Trans.Elem1, eip1);
2321  }
2322  a *= rho_p;
2323  b *= rho_p;
2324  }
2325 
2326  w = ip.weight * (a+b);
2327  if (w != 0.0)
2328  {
2329  for (int i = 0; i < ndof1; i++)
2330  for (int j = 0; j < ndof1; j++)
2331  {
2332  elmat(i, j) += w * shape1(i) * shape1(j);
2333  }
2334  }
2335 
2336  if (ndof2)
2337  {
2338  el2.CalcShape(eip2, shape2);
2339 
2340  if (w != 0.0)
2341  for (int i = 0; i < ndof2; i++)
2342  for (int j = 0; j < ndof1; j++)
2343  {
2344  elmat(ndof1+i, j) -= w * shape2(i) * shape1(j);
2345  }
2346 
2347  w = ip.weight * (b-a);
2348  if (w != 0.0)
2349  {
2350  for (int i = 0; i < ndof2; i++)
2351  for (int j = 0; j < ndof2; j++)
2352  {
2353  elmat(ndof1+i, ndof1+j) += w * shape2(i) * shape2(j);
2354  }
2355 
2356  for (int i = 0; i < ndof1; i++)
2357  for (int j = 0; j < ndof2; j++)
2358  {
2359  elmat(i, ndof1+j) -= w * shape1(i) * shape2(j);
2360  }
2361  }
2362  }
2363  }
2364 }
2365 
2366 void DGDiffusionIntegrator::AssembleFaceMatrix(
2367  const FiniteElement &el1, const FiniteElement &el2,
2368  FaceElementTransformations &Trans, DenseMatrix &elmat)
2369 {
2370  int dim, ndof1, ndof2, ndofs;
2371  bool kappa_is_nonzero = (kappa != 0.);
2372  double w, wq = 0.0;
2373 
2374  dim = el1.GetDim();
2375  ndof1 = el1.GetDof();
2376 
2377  nor.SetSize(dim);
2378  nh.SetSize(dim);
2379  ni.SetSize(dim);
2380  adjJ.SetSize(dim);
2381  if (MQ)
2382  {
2383  mq.SetSize(dim);
2384  }
2385 
2386  shape1.SetSize(ndof1);
2387  dshape1.SetSize(ndof1, dim);
2388  dshape1dn.SetSize(ndof1);
2389  if (Trans.Elem2No >= 0)
2390  {
2391  ndof2 = el2.GetDof();
2392  shape2.SetSize(ndof2);
2393  dshape2.SetSize(ndof2, dim);
2394  dshape2dn.SetSize(ndof2);
2395  }
2396  else
2397  {
2398  ndof2 = 0;
2399  }
2400 
2401  ndofs = ndof1 + ndof2;
2402  elmat.SetSize(ndofs);
2403  elmat = 0.0;
2404  if (kappa_is_nonzero)
2405  {
2406  jmat.SetSize(ndofs);
2407  jmat = 0.;
2408  }
2409 
2410  const IntegrationRule *ir = IntRule;
2411  if (ir == NULL)
2412  {
2413  // a simple choice for the integration order; is this OK?
2414  int order;
2415  if (ndof2)
2416  {
2417  order = 2*max(el1.GetOrder(), el2.GetOrder());
2418  }
2419  else
2420  {
2421  order = 2*el1.GetOrder();
2422  }
2423  ir = &IntRules.Get(Trans.FaceGeom, order);
2424  }
2425 
2426  // assemble: < {(Q \nabla u).n},[v] > --> elmat
2427  // kappa < {h^{-1} Q} [u],[v] > --> jmat
2428  for (int p = 0; p < ir->GetNPoints(); p++)
2429  {
2430  const IntegrationPoint &ip = ir->IntPoint(p);
2431  IntegrationPoint eip1, eip2;
2432 
2433  Trans.Loc1.Transform(ip, eip1);
2434  Trans.Face->SetIntPoint(&ip);
2435  if (dim == 1)
2436  {
2437  nor(0) = 2*eip1.x - 1.0;
2438  }
2439  else
2440  {
2441  CalcOrtho(Trans.Face->Jacobian(), nor);
2442  }
2443 
2444  el1.CalcShape(eip1, shape1);
2445  el1.CalcDShape(eip1, dshape1);
2446  Trans.Elem1->SetIntPoint(&eip1);
2447  w = ip.weight/Trans.Elem1->Weight();
2448  if (ndof2)
2449  {
2450  w /= 2;
2451  }
2452  if (!MQ)
2453  {
2454  if (Q)
2455  {
2456  w *= Q->Eval(*Trans.Elem1, eip1);
2457  }
2458  ni.Set(w, nor);
2459  }
2460  else
2461  {
2462  nh.Set(w, nor);
2463  MQ->Eval(mq, *Trans.Elem1, eip1);
2464  mq.MultTranspose(nh, ni);
2465  }
2466  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
2467  adjJ.Mult(ni, nh);
2468  if (kappa_is_nonzero)
2469  {
2470  wq = ni * nor;
2471  }
2472  // Note: in the jump term, we use 1/h1 = |nor|/det(J1) which is
2473  // independent of Loc1 and always gives the size of element 1 in
2474  // direction perpendicular to the face. Indeed, for linear transformation
2475  // |nor|=measure(face)/measure(ref. face),
2476  // det(J1)=measure(element)/measure(ref. element),
2477  // and the ratios measure(ref. element)/measure(ref. face) are
2478  // compatible for all element/face pairs.
2479  // For example: meas(ref. tetrahedron)/meas(ref. triangle) = 1/3, and
2480  // for any tetrahedron vol(tet)=(1/3)*height*area(base).
2481  // For interior faces: q_e/h_e=(q1/h1+q2/h2)/2.
2482 
2483  dshape1.Mult(nh, dshape1dn);
2484  for (int i = 0; i < ndof1; i++)
2485  for (int j = 0; j < ndof1; j++)
2486  {
2487  elmat(i, j) += shape1(i) * dshape1dn(j);
2488  }
2489 
2490  if (ndof2)
2491  {
2492  Trans.Loc2.Transform(ip, eip2);
2493  el2.CalcShape(eip2, shape2);
2494  el2.CalcDShape(eip2, dshape2);
2495  Trans.Elem2->SetIntPoint(&eip2);
2496  w = ip.weight/2/Trans.Elem2->Weight();
2497  if (!MQ)
2498  {
2499  if (Q)
2500  {
2501  w *= Q->Eval(*Trans.Elem2, eip2);
2502  }
2503  ni.Set(w, nor);
2504  }
2505  else
2506  {
2507  nh.Set(w, nor);
2508  MQ->Eval(mq, *Trans.Elem2, eip2);
2509  mq.MultTranspose(nh, ni);
2510  }
2511  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
2512  adjJ.Mult(ni, nh);
2513  if (kappa_is_nonzero)
2514  {
2515  wq += ni * nor;
2516  }
2517 
2518  dshape2.Mult(nh, dshape2dn);
2519 
2520  for (int i = 0; i < ndof1; i++)
2521  for (int j = 0; j < ndof2; j++)
2522  {
2523  elmat(i, ndof1 + j) += shape1(i) * dshape2dn(j);
2524  }
2525 
2526  for (int i = 0; i < ndof2; i++)
2527  for (int j = 0; j < ndof1; j++)
2528  {
2529  elmat(ndof1 + i, j) -= shape2(i) * dshape1dn(j);
2530  }
2531 
2532  for (int i = 0; i < ndof2; i++)
2533  for (int j = 0; j < ndof2; j++)
2534  {
2535  elmat(ndof1 + i, ndof1 + j) -= shape2(i) * dshape2dn(j);
2536  }
2537  }
2538 
2539  if (kappa_is_nonzero)
2540  {
2541  // only assemble the lower triangular part of jmat
2542  wq *= kappa;
2543  for (int i = 0; i < ndof1; i++)
2544  {
2545  const double wsi = wq*shape1(i);
2546  for (int j = 0; j <= i; j++)
2547  {
2548  jmat(i, j) += wsi * shape1(j);
2549  }
2550  }
2551  if (ndof2)
2552  {
2553  for (int i = 0; i < ndof2; i++)
2554  {
2555  const int i2 = ndof1 + i;
2556  const double wsi = wq*shape2(i);
2557  for (int j = 0; j < ndof1; j++)
2558  {
2559  jmat(i2, j) -= wsi * shape1(j);
2560  }
2561  for (int j = 0; j <= i; j++)
2562  {
2563  jmat(i2, ndof1 + j) += wsi * shape2(j);
2564  }
2565  }
2566  }
2567  }
2568  }
2569 
2570  // elmat := -elmat + sigma*elmat^t + jmat
2571  if (kappa_is_nonzero)
2572  {
2573  for (int i = 0; i < ndofs; i++)
2574  {
2575  for (int j = 0; j < i; j++)
2576  {
2577  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
2578  elmat(i,j) = sigma*aji - aij + mij;
2579  elmat(j,i) = sigma*aij - aji + mij;
2580  }
2581  elmat(i,i) = (sigma - 1.)*elmat(i,i) + jmat(i,i);
2582  }
2583  }
2584  else
2585  {
2586  for (int i = 0; i < ndofs; i++)
2587  {
2588  for (int j = 0; j < i; j++)
2589  {
2590  double aij = elmat(i,j), aji = elmat(j,i);
2591  elmat(i,j) = sigma*aji - aij;
2592  elmat(j,i) = sigma*aij - aji;
2593  }
2594  elmat(i,i) *= (sigma - 1.);
2595  }
2596  }
2597 }
2598 
2599 
2600 // static method
2601 void DGElasticityIntegrator::AssembleBlock(
2602  const int dim, const int row_ndofs, const int col_ndofs,
2603  const int row_offset, const int col_offset,
2604  const double jmatcoef, const Vector &col_nL, const Vector &col_nM,
2605  const Vector &row_shape, const Vector &col_shape,
2606  const Vector &col_dshape_dnM, const DenseMatrix &col_dshape,
2607  DenseMatrix &elmat, DenseMatrix &jmat)
2608 {
2609  for (int jm = 0, j = col_offset; jm < dim; ++jm)
2610  {
2611  for (int jdof = 0; jdof < col_ndofs; ++jdof, ++j)
2612  {
2613  const double t2 = col_dshape_dnM(jdof);
2614  for (int im = 0, i = row_offset; im < dim; ++im)
2615  {
2616  const double t1 = col_dshape(jdof, jm) * col_nL(im);
2617  const double t3 = col_dshape(jdof, im) * col_nM(jm);
2618  const double tt = t1 + ((im == jm) ? t2 : 0.0) + t3;
2619  for (int idof = 0; idof < row_ndofs; ++idof, ++i)
2620  {
2621  elmat(i, j) += row_shape(idof) * tt;
2622  }
2623  }
2624  }
2625  }
2626 
2627  if (jmatcoef == 0.0) { return; }
2628 
2629  for (int d = 0; d < dim; ++d)
2630  {
2631  const int jo = col_offset + d*col_ndofs;
2632  const int io = row_offset + d*row_ndofs;
2633  for (int jdof = 0, j = jo; jdof < col_ndofs; ++jdof, ++j)
2634  {
2635  const double sj = jmatcoef * col_shape(jdof);
2636  for (int i = max(io,j), idof = i - io; idof < row_ndofs; ++idof, ++i)
2637  {
2638  jmat(i, j) += row_shape(idof) * sj;
2639  }
2640  }
2641  }
2642 }
2643 
2644 void DGElasticityIntegrator::AssembleFaceMatrix(
2645  const FiniteElement &el1, const FiniteElement &el2,
2646  FaceElementTransformations &Trans, DenseMatrix &elmat)
2647 {
2648 #ifdef MFEM_THREAD_SAFE
2649  // For descriptions of these variables, see the class declaration.
2650  Vector shape1, shape2;
2651  DenseMatrix dshape1, dshape2;
2652  DenseMatrix adjJ;
2653  DenseMatrix dshape1_ps, dshape2_ps;
2654  Vector nor;
2655  Vector nL1, nL2;
2656  Vector nM1, nM2;
2657  Vector dshape1_dnM, dshape2_dnM;
2658  DenseMatrix jmat;
2659 #endif
2660 
2661  const int dim = el1.GetDim();
2662  const int ndofs1 = el1.GetDof();
2663  const int ndofs2 = (Trans.Elem2No >= 0) ? el2.GetDof() : 0;
2664  const int nvdofs = dim*(ndofs1 + ndofs2);
2665 
2666  // Initially 'elmat' corresponds to the term:
2667  // < { sigma(u) . n }, [v] > =
2668  // < { (lambda div(u) I + mu (grad(u) + grad(u)^T)) . n }, [v] >
2669  // But eventually, it's going to be replaced by:
2670  // elmat := -elmat + alpha*elmat^T + jmat
2671  elmat.SetSize(nvdofs);
2672  elmat = 0.;
2673 
2674  const bool kappa_is_nonzero = (kappa != 0.0);
2675  if (kappa_is_nonzero)
2676  {
2677  jmat.SetSize(nvdofs);
2678  jmat = 0.;
2679  }
2680 
2681  adjJ.SetSize(dim);
2682  shape1.SetSize(ndofs1);
2683  dshape1.SetSize(ndofs1, dim);
2684  dshape1_ps.SetSize(ndofs1, dim);
2685  nor.SetSize(dim);
2686  nL1.SetSize(dim);
2687  nM1.SetSize(dim);
2688  dshape1_dnM.SetSize(ndofs1);
2689 
2690  if (ndofs2)
2691  {
2692  shape2.SetSize(ndofs2);
2693  dshape2.SetSize(ndofs2, dim);
2694  dshape2_ps.SetSize(ndofs2, dim);
2695  nL2.SetSize(dim);
2696  nM2.SetSize(dim);
2697  dshape2_dnM.SetSize(ndofs2);
2698  }
2699 
2700  const IntegrationRule *ir = IntRule;
2701  if (ir == NULL)
2702  {
2703  // a simple choice for the integration order; is this OK?
2704  const int order = 2 * max(el1.GetOrder(), ndofs2 ? el2.GetOrder() : 0);
2705  ir = &IntRules.Get(Trans.FaceGeom, order);
2706  }
2707 
2708  for (int pind = 0; pind < ir->GetNPoints(); ++pind)
2709  {
2710  const IntegrationPoint &ip = ir->IntPoint(pind);
2711  IntegrationPoint eip1, eip2; // integration point in the reference space
2712  Trans.Loc1.Transform(ip, eip1);
2713  Trans.Face->SetIntPoint(&ip);
2714  Trans.Elem1->SetIntPoint(&eip1);
2715 
2716  el1.CalcShape(eip1, shape1);
2717  el1.CalcDShape(eip1, dshape1);
2718 
2719  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
2720  Mult(dshape1, adjJ, dshape1_ps);
2721 
2722  if (dim == 1)
2723  {
2724  nor(0) = 2*eip1.x - 1.0;
2725  }
2726  else
2727  {
2728  CalcOrtho(Trans.Face->Jacobian(), nor);
2729  }
2730 
2731  double w, wLM;
2732  if (ndofs2)
2733  {
2734  Trans.Loc2.Transform(ip, eip2);
2735  Trans.Elem2->SetIntPoint(&eip2);
2736  el2.CalcShape(eip2, shape2);
2737  el2.CalcDShape(eip2, dshape2);
2738  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
2739  Mult(dshape2, adjJ, dshape2_ps);
2740 
2741  w = ip.weight/2;
2742  const double w2 = w / Trans.Elem2->Weight();
2743  const double wL2 = w2 * lambda->Eval(*Trans.Elem2, eip2);
2744  const double wM2 = w2 * mu->Eval(*Trans.Elem2, eip2);
2745  nL2.Set(wL2, nor);
2746  nM2.Set(wM2, nor);
2747  wLM = (wL2 + 2.0*wM2);
2748  dshape2_ps.Mult(nM2, dshape2_dnM);
2749  }
2750  else
2751  {
2752  w = ip.weight;
2753  wLM = 0.0;
2754  }
2755 
2756  {
2757  const double w1 = w / Trans.Elem1->Weight();
2758  const double wL1 = w1 * lambda->Eval(*Trans.Elem1, eip1);
2759  const double wM1 = w1 * mu->Eval(*Trans.Elem1, eip1);
2760  nL1.Set(wL1, nor);
2761  nM1.Set(wM1, nor);
2762  wLM += (wL1 + 2.0*wM1);
2763  dshape1_ps.Mult(nM1, dshape1_dnM);
2764  }
2765 
2766  const double jmatcoef = kappa * (nor*nor) * wLM;
2767 
2768  // (1,1) block
2769  AssembleBlock(
2770  dim, ndofs1, ndofs1, 0, 0, jmatcoef, nL1, nM1,
2771  shape1, shape1, dshape1_dnM, dshape1_ps, elmat, jmat);
2772 
2773  if (ndofs2 == 0) { continue; }
2774 
2775  // In both elmat and jmat, shape2 appears only with a minus sign.
2776  shape2.Neg();
2777 
2778  // (1,2) block
2779  AssembleBlock(
2780  dim, ndofs1, ndofs2, 0, dim*ndofs1, jmatcoef, nL2, nM2,
2781  shape1, shape2, dshape2_dnM, dshape2_ps, elmat, jmat);
2782  // (2,1) block
2783  AssembleBlock(
2784  dim, ndofs2, ndofs1, dim*ndofs1, 0, jmatcoef, nL1, nM1,
2785  shape2, shape1, dshape1_dnM, dshape1_ps, elmat, jmat);
2786  // (2,2) block
2787  AssembleBlock(
2788  dim, ndofs2, ndofs2, dim*ndofs1, dim*ndofs1, jmatcoef, nL2, nM2,
2789  shape2, shape2, dshape2_dnM, dshape2_ps, elmat, jmat);
2790  }
2791 
2792  // elmat := -elmat + alpha*elmat^t + jmat
2793  if (kappa_is_nonzero)
2794  {
2795  for (int i = 0; i < nvdofs; ++i)
2796  {
2797  for (int j = 0; j < i; ++j)
2798  {
2799  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
2800  elmat(i,j) = alpha*aji - aij + mij;
2801  elmat(j,i) = alpha*aij - aji + mij;
2802  }
2803  elmat(i,i) = (alpha - 1.)*elmat(i,i) + jmat(i,i);
2804  }
2805  }
2806  else
2807  {
2808  for (int i = 0; i < nvdofs; ++i)
2809  {
2810  for (int j = 0; j < i; ++j)
2811  {
2812  double aij = elmat(i,j), aji = elmat(j,i);
2813  elmat(i,j) = alpha*aji - aij;
2814  elmat(j,i) = alpha*aij - aji;
2815  }
2816  elmat(i,i) *= (alpha - 1.);
2817  }
2818  }
2819 }
2820 
2821 
2822 void TraceJumpIntegrator::AssembleFaceMatrix(
2823  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2824  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2825  DenseMatrix &elmat)
2826 {
2827  int i, j, face_ndof, ndof1, ndof2;
2828  int order;
2829 
2830  double w;
2831 
2832  face_ndof = trial_face_fe.GetDof();
2833  ndof1 = test_fe1.GetDof();
2834 
2835  face_shape.SetSize(face_ndof);
2836  shape1.SetSize(ndof1);
2837 
2838  if (Trans.Elem2No >= 0)
2839  {
2840  ndof2 = test_fe2.GetDof();
2841  shape2.SetSize(ndof2);
2842  }
2843  else
2844  {
2845  ndof2 = 0;
2846  }
2847 
2848  elmat.SetSize(ndof1 + ndof2, face_ndof);
2849  elmat = 0.0;
2850 
2851  const IntegrationRule *ir = IntRule;
2852  if (ir == NULL)
2853  {
2854  if (Trans.Elem2No >= 0)
2855  {
2856  order = max(test_fe1.GetOrder(), test_fe2.GetOrder());
2857  }
2858  else
2859  {
2860  order = test_fe1.GetOrder();
2861  }
2862  order += trial_face_fe.GetOrder();
2863  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2864  {
2865  order += Trans.Face->OrderW();
2866  }
2867  ir = &IntRules.Get(Trans.FaceGeom, order);
2868  }
2869 
2870  for (int p = 0; p < ir->GetNPoints(); p++)
2871  {
2872  const IntegrationPoint &ip = ir->IntPoint(p);
2873  IntegrationPoint eip1, eip2;
2874  // Trace finite element shape function
2875  Trans.Face->SetIntPoint(&ip);
2876  trial_face_fe.CalcShape(ip, face_shape);
2877  // Side 1 finite element shape function
2878  Trans.Loc1.Transform(ip, eip1);
2879  test_fe1.CalcShape(eip1, shape1);
2880  Trans.Elem1->SetIntPoint(&eip1);
2881  if (ndof2)
2882  {
2883  // Side 2 finite element shape function
2884  Trans.Loc2.Transform(ip, eip2);
2885  test_fe2.CalcShape(eip2, shape2);
2886  Trans.Elem2->SetIntPoint(&eip2);
2887  }
2888  w = ip.weight;
2889  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2890  {
2891  w *= Trans.Face->Weight();
2892  }
2893  face_shape *= w;
2894  for (i = 0; i < ndof1; i++)
2895  for (j = 0; j < face_ndof; j++)
2896  {
2897  elmat(i, j) += shape1(i) * face_shape(j);
2898  }
2899  if (ndof2)
2900  {
2901  // Subtract contribution from side 2
2902  for (i = 0; i < ndof2; i++)
2903  for (j = 0; j < face_ndof; j++)
2904  {
2905  elmat(ndof1+i, j) -= shape2(i) * face_shape(j);
2906  }
2907  }
2908  }
2909 }
2910 
2911 void NormalTraceJumpIntegrator::AssembleFaceMatrix(
2912  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2913  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2914  DenseMatrix &elmat)
2915 {
2916  int i, j, face_ndof, ndof1, ndof2, dim;
2917  int order;
2918 
2919  MFEM_VERIFY(trial_face_fe.GetMapType() == FiniteElement::VALUE, "");
2920 
2921  face_ndof = trial_face_fe.GetDof();
2922  ndof1 = test_fe1.GetDof();
2923  dim = test_fe1.GetDim();
2924 
2925  face_shape.SetSize(face_ndof);
2926  normal.SetSize(dim);
2927  shape1.SetSize(ndof1,dim);
2928  shape1_n.SetSize(ndof1);
2929 
2930  if (Trans.Elem2No >= 0)
2931  {
2932  ndof2 = test_fe2.GetDof();
2933  shape2.SetSize(ndof2,dim);
2934  shape2_n.SetSize(ndof2);
2935  }
2936  else
2937  {
2938  ndof2 = 0;
2939  }
2940 
2941  elmat.SetSize(ndof1 + ndof2, face_ndof);
2942  elmat = 0.0;
2943 
2944  const IntegrationRule *ir = IntRule;
2945  if (ir == NULL)
2946  {
2947  if (Trans.Elem2No >= 0)
2948  {
2949  order = max(test_fe1.GetOrder(), test_fe2.GetOrder()) - 1;
2950  }
2951  else
2952  {
2953  order = test_fe1.GetOrder() - 1;
2954  }
2955  order += trial_face_fe.GetOrder();
2956  ir = &IntRules.Get(Trans.FaceGeom, order);
2957  }
2958 
2959  for (int p = 0; p < ir->GetNPoints(); p++)
2960  {
2961  const IntegrationPoint &ip = ir->IntPoint(p);
2962  IntegrationPoint eip1, eip2;
2963  // Trace finite element shape function
2964  trial_face_fe.CalcShape(ip, face_shape);
2965  Trans.Loc1.Transf.SetIntPoint(&ip);
2966  CalcOrtho(Trans.Loc1.Transf.Jacobian(), normal);
2967  // Side 1 finite element shape function
2968  Trans.Loc1.Transform(ip, eip1);
2969  test_fe1.CalcVShape(eip1, shape1);
2970  shape1.Mult(normal, shape1_n);
2971  if (ndof2)
2972  {
2973  // Side 2 finite element shape function
2974  Trans.Loc2.Transform(ip, eip2);
2975  test_fe2.CalcVShape(eip2, shape2);
2976  Trans.Loc2.Transf.SetIntPoint(&ip);
2977  CalcOrtho(Trans.Loc2.Transf.Jacobian(), normal);
2978  shape2.Mult(normal, shape2_n);
2979  }
2980  face_shape *= ip.weight;
2981  for (i = 0; i < ndof1; i++)
2982  for (j = 0; j < face_ndof; j++)
2983  {
2984  elmat(i, j) -= shape1_n(i) * face_shape(j);
2985  }
2986  if (ndof2)
2987  {
2988  // Subtract contribution from side 2
2989  for (i = 0; i < ndof2; i++)
2990  for (j = 0; j < face_ndof; j++)
2991  {
2992  elmat(ndof1+i, j) += shape2_n(i) * face_shape(j);
2993  }
2994  }
2995  }
2996 }
2997 
2998 
2999 void NormalInterpolator::AssembleElementMatrix2(
3000  const FiniteElement &dom_fe, const FiniteElement &ran_fe,
3001  ElementTransformation &Trans, DenseMatrix &elmat)
3002 {
3003  int spaceDim = Trans.GetSpaceDim();
3004  elmat.SetSize(ran_fe.GetDof(), spaceDim*dom_fe.GetDof());
3005  Vector n(spaceDim), shape(dom_fe.GetDof());
3006 
3007  const IntegrationRule &ran_nodes = ran_fe.GetNodes();
3008  for (int i = 0; i < ran_nodes.Size(); i++)
3009  {
3010  const IntegrationPoint &ip = ran_nodes.IntPoint(i);
3011  Trans.SetIntPoint(&ip);
3012  CalcOrtho(Trans.Jacobian(), n);
3013  dom_fe.CalcShape(ip, shape);
3014  for (int j = 0; j < shape.Size(); j++)
3015  {
3016  for (int d = 0; d < spaceDim; d++)
3017  {
3018  elmat(i, j+d*shape.Size()) = shape(j)*n(d);
3019  }
3020  }
3021  }
3022 }
3023 
3024 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:222
Abstract class for Finite Elements.
Definition: fe.hpp:46
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:3216
const DenseMatrix & AdjugateJacobian()
Definition: eltrans.hpp:69
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:115
ElementTransformation * Face
Definition: eltrans.hpp:158
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3618
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:94
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe.cpp:38
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:3596
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:861
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:310
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:271
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:51
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2902
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:124
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:159
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:127
void AddMult_a_ABt(double a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
Definition: densemat.cpp:3432
int GetMapType() const
Definition: fe.hpp:133
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:90
ElementTransformation * Elem2
Definition: eltrans.hpp:158
double * GetData() const
Definition: vector.hpp:114
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3132
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:148
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:225
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:193
int dim
Definition: ex3.cpp:47
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:67
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:159
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:3639
const DenseMatrix & Jacobian()
Definition: eltrans.hpp:64
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
For scalar fields; preserves point values.
Definition: fe.hpp:84
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:566
const IntegrationRule & GetNodes() const
Definition: fe.hpp:166
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:3661
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:3010
int GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:118
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:2292
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:3587
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:3339
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:121
void mfem_error(const char *msg)
Definition: error.cpp:106
void ClearExternalData()
Definition: densemat.hpp:72
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:256
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:3158
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe.cpp:52
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0&lt;=i&lt;A.Height, 0&lt;=j&lt;A.Width.
Definition: densemat.cpp:2557
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:87
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:223
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
Definition: fe.cpp:66
void AddMultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
Definition: densemat.cpp:3396
Class for integration point with weight.
Definition: intrules.hpp:25
virtual int OrderGrad(const FiniteElement *fe)=0
order of adj(J)^t.grad(fi)
IsoparametricTransformation Transf
Definition: eltrans.hpp:149
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:2356
ElementTransformation * Elem1
Definition: eltrans.hpp:158
const double alpha
Definition: ex15.cpp:337
double kappa
Definition: ex3.cpp:46
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:3491
Vector data type.
Definition: vector.hpp:36
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:346
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
Definition: densemat.cpp:3172
virtual int GetSpaceDim()=0
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:82
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
int GetRangeType() const
Definition: fe.hpp:129
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:343
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:3548
void Neg()
(*this) = -(*this)
Definition: vector.cpp:256
const int ir_order