MFEM  v4.4.0
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-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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::AssemblePA(const FiniteElementSpace&)
24 {
25  mfem_error ("BilinearFormIntegrator::AssemblePA(fes)\n"
26  " is not implemented for this class.");
27 }
28 
29 void BilinearFormIntegrator::AssemblePA(const FiniteElementSpace&,
30  const FiniteElementSpace&)
31 {
32  mfem_error ("BilinearFormIntegrator::AssemblePA(fes, fes)\n"
33  " is not implemented for this class.");
34 }
35 
36 void BilinearFormIntegrator::AssemblePAInteriorFaces(const FiniteElementSpace&)
37 {
38  mfem_error ("BilinearFormIntegrator::AssemblePAInteriorFaces(...)\n"
39  " is not implemented for this class.");
40 }
41 
42 void BilinearFormIntegrator::AssemblePABoundaryFaces(const FiniteElementSpace&)
43 {
44  mfem_error ("BilinearFormIntegrator::AssemblePABoundaryFaces(...)\n"
45  " is not implemented for this class.");
46 }
47 
48 void BilinearFormIntegrator::AssembleDiagonalPA(Vector &)
49 {
50  mfem_error ("BilinearFormIntegrator::AssembleDiagonalPA(...)\n"
51  " is not implemented for this class.");
52 }
53 
54 void BilinearFormIntegrator::AssembleEA(const FiniteElementSpace &fes,
55  Vector &emat,
56  const bool add)
57 {
58  mfem_error ("BilinearFormIntegrator::AssembleEA(...)\n"
59  " is not implemented for this class.");
60 }
61 
62 void BilinearFormIntegrator::AssembleEAInteriorFaces(const FiniteElementSpace
63  &fes,
64  Vector &ea_data_int,
65  Vector &ea_data_ext,
66  const bool add)
67 {
68  mfem_error ("BilinearFormIntegrator::AssembleEAInteriorFaces(...)\n"
69  " is not implemented for this class.");
70 }
71 
72 void BilinearFormIntegrator::AssembleEABoundaryFaces(const FiniteElementSpace
73  &fes,
74  Vector &ea_data_bdr,
75  const bool add)
76 {
77  mfem_error ("BilinearFormIntegrator::AssembleEABoundaryFaces(...)\n"
78  " is not implemented for this class.");
79 }
80 
81 void BilinearFormIntegrator::AssembleDiagonalPA_ADAt(const Vector &, Vector &)
82 {
83  MFEM_ABORT("BilinearFormIntegrator::AssembleDiagonalPA_ADAt(...)\n"
84  " is not implemented for this class.");
85 }
86 
87 void BilinearFormIntegrator::AddMultPA(const Vector &, Vector &) const
88 {
89  mfem_error ("BilinearFormIntegrator::MultAssembled(...)\n"
90  " is not implemented for this class.");
91 }
92 
93 void BilinearFormIntegrator::AddMultTransposePA(const Vector &, Vector &) const
94 {
95  mfem_error ("BilinearFormIntegrator::AddMultTransposePA(...)\n"
96  " is not implemented for this class.");
97 }
98 
99 void BilinearFormIntegrator::AssembleMF(const FiniteElementSpace &fes)
100 {
101  mfem_error ("BilinearFormIntegrator::AssembleMF(...)\n"
102  " is not implemented for this class.");
103 }
104 
105 void BilinearFormIntegrator::AddMultMF(const Vector &, Vector &) const
106 {
107  mfem_error ("BilinearFormIntegrator::AddMultMF(...)\n"
108  " is not implemented for this class.");
109 }
110 
111 void BilinearFormIntegrator::AddMultTransposeMF(const Vector &, Vector &) const
112 {
113  mfem_error ("BilinearFormIntegrator::AddMultTransposeMF(...)\n"
114  " is not implemented for this class.");
115 }
116 
117 void BilinearFormIntegrator::AssembleDiagonalMF(Vector &)
118 {
119  mfem_error ("BilinearFormIntegrator::AssembleDiagonalMF(...)\n"
120  " is not implemented for this class.");
121 }
122 
123 void BilinearFormIntegrator::AssembleElementMatrix (
125  DenseMatrix &elmat )
126 {
127  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix(...)\n"
128  " is not implemented for this class.");
129 }
130 
131 void BilinearFormIntegrator::AssembleElementMatrix2 (
132  const FiniteElement &el1, const FiniteElement &el2,
134 {
135  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix2(...)\n"
136  " is not implemented for this class.");
137 }
138 
139 void BilinearFormIntegrator::AssembleFaceMatrix (
140  const FiniteElement &el1, const FiniteElement &el2,
142 {
143  mfem_error ("BilinearFormIntegrator::AssembleFaceMatrix(...)\n"
144  " is not implemented for this class.");
145 }
146 
147 void BilinearFormIntegrator::AssembleFaceMatrix(
148  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
150  DenseMatrix &elmat)
151 {
152  MFEM_ABORT("AssembleFaceMatrix (mixed form) is not implemented for this"
153  " Integrator class.");
154 }
155 
156 void BilinearFormIntegrator::AssembleElementVector(
157  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
158  Vector &elvect)
159 {
160  // Note: This default implementation is general but not efficient
161  DenseMatrix elmat;
162  AssembleElementMatrix(el, Tr, elmat);
163  elvect.SetSize(elmat.Height());
164  elmat.Mult(elfun, elvect);
165 }
166 
167 void BilinearFormIntegrator::AssembleFaceVector(
168  const FiniteElement &el1, const FiniteElement &el2,
169  FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
170 {
171  // Note: This default implementation is general but not efficient
172  DenseMatrix elmat;
173  AssembleFaceMatrix(el1, el2, Tr, elmat);
174  elvect.SetSize(elmat.Height());
175  elmat.Mult(elfun, elvect);
176 }
177 
178 void TransposeIntegrator::SetIntRule(const IntegrationRule *ir)
179 {
180  IntRule = ir;
181  bfi->SetIntRule(ir);
182 }
183 
184 void TransposeIntegrator::AssembleElementMatrix (
186 {
187  bfi -> AssembleElementMatrix (el, Trans, bfi_elmat);
188  // elmat = bfi_elmat^t
189  elmat.Transpose (bfi_elmat);
190 }
191 
192 void TransposeIntegrator::AssembleElementMatrix2 (
193  const FiniteElement &trial_fe, const FiniteElement &test_fe,
195 {
196  bfi -> AssembleElementMatrix2 (test_fe, trial_fe, Trans, bfi_elmat);
197  // elmat = bfi_elmat^t
198  elmat.Transpose (bfi_elmat);
199 }
200 
201 void TransposeIntegrator::AssembleFaceMatrix (
202  const FiniteElement &el1, const FiniteElement &el2,
204 {
205  bfi -> AssembleFaceMatrix (el1, el2, Trans, bfi_elmat);
206  // elmat = bfi_elmat^t
207  elmat.Transpose (bfi_elmat);
208 }
209 
210 void LumpedIntegrator::SetIntRule(const IntegrationRule *ir)
211 {
212  IntRule = ir;
213  bfi->SetIntRule(ir);
214 }
215 
216 void LumpedIntegrator::AssembleElementMatrix (
218 {
219  bfi -> AssembleElementMatrix (el, Trans, elmat);
220  elmat.Lump();
221 }
222 
223 void InverseIntegrator::SetIntRule(const IntegrationRule *ir)
224 {
225  IntRule = ir;
226  integrator->SetIntRule(ir);
227 }
228 
229 void InverseIntegrator::AssembleElementMatrix(
231 {
232  integrator->AssembleElementMatrix(el, Trans, elmat);
233  elmat.Invert();
234 }
235 
236 void SumIntegrator::SetIntRule(const IntegrationRule *ir)
237 {
238  IntRule = ir;
239  for (int i = 0; i < integrators.Size(); i++)
240  {
241  integrators[i]->SetIntRule(ir);
242  }
243 }
244 
245 void SumIntegrator::AssembleElementMatrix(
247 {
248  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
249 
250  integrators[0]->AssembleElementMatrix(el, Trans, elmat);
251  for (int i = 1; i < integrators.Size(); i++)
252  {
253  integrators[i]->AssembleElementMatrix(el, Trans, elem_mat);
254  elmat += elem_mat;
255  }
256 }
257 
258 void SumIntegrator::AssembleElementMatrix2(
259  const FiniteElement &el1, const FiniteElement &el2,
261 {
262  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
263 
264  integrators[0]->AssembleElementMatrix2(el1, el2, Trans, elmat);
265  for (int i = 1; i < integrators.Size(); i++)
266  {
267  integrators[i]->AssembleElementMatrix2(el1, el2, Trans, elem_mat);
268  elmat += elem_mat;
269  }
270 }
271 
272 void SumIntegrator::AssembleFaceMatrix(
273  const FiniteElement &el1, const FiniteElement &el2,
275 {
276  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
277 
278  integrators[0]->AssembleFaceMatrix(el1, el2, Trans, elmat);
279  for (int i = 1; i < integrators.Size(); i++)
280  {
281  integrators[i]->AssembleFaceMatrix(el1, el2, Trans, elem_mat);
282  elmat += elem_mat;
283  }
284 }
285 
286 void SumIntegrator::AssembleFaceMatrix(
287  const FiniteElement &tr_fe,
288  const FiniteElement &te_fe1, const FiniteElement &te_fe2,
290 {
291  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
292 
293  integrators[0]->AssembleFaceMatrix(tr_fe, te_fe1, te_fe2, Trans, elmat);
294  for (int i = 1; i < integrators.Size(); i++)
295  {
296  integrators[i]->AssembleFaceMatrix(tr_fe, te_fe1, te_fe2, Trans, elem_mat);
297  elmat += elem_mat;
298  }
299 }
300 
301 void SumIntegrator::AssemblePA(const FiniteElementSpace& fes)
302 {
303  for (int i = 0; i < integrators.Size(); i++)
304  {
305  integrators[i]->AssemblePA(fes);
306  }
307 }
308 
309 void SumIntegrator::AssembleDiagonalPA(Vector &diag)
310 {
311  for (int i = 0; i < integrators.Size(); i++)
312  {
313  integrators[i]->AssembleDiagonalPA(diag);
314  }
315 }
316 
317 void SumIntegrator::AssemblePAInteriorFaces(const FiniteElementSpace &fes)
318 {
319  for (int i = 0; i < integrators.Size(); i++)
320  {
321  integrators[i]->AssemblePAInteriorFaces(fes);
322  }
323 }
324 
325 void SumIntegrator::AssemblePABoundaryFaces(const FiniteElementSpace &fes)
326 {
327  for (int i = 0; i < integrators.Size(); i++)
328  {
329  integrators[i]->AssemblePABoundaryFaces(fes);
330  }
331 }
332 
333 void SumIntegrator::AddMultPA(const Vector& x, Vector& y) const
334 {
335  for (int i = 0; i < integrators.Size(); i++)
336  {
337  integrators[i]->AddMultPA(x, y);
338  }
339 }
340 
341 void SumIntegrator::AddMultTransposePA(const Vector &x, Vector &y) const
342 {
343  for (int i = 0; i < integrators.Size(); i++)
344  {
345  integrators[i]->AddMultTransposePA(x, y);
346  }
347 }
348 
349 void SumIntegrator::AssembleMF(const FiniteElementSpace &fes)
350 {
351  for (int i = 0; i < integrators.Size(); i++)
352  {
353  integrators[i]->AssembleMF(fes);
354  }
355 }
356 
357 void SumIntegrator::AddMultMF(const Vector& x, Vector& y) const
358 {
359  for (int i = 0; i < integrators.Size(); i++)
360  {
361  integrators[i]->AddMultTransposeMF(x, y);
362  }
363 }
364 
365 void SumIntegrator::AddMultTransposeMF(const Vector &x, Vector &y) const
366 {
367  for (int i = 0; i < integrators.Size(); i++)
368  {
369  integrators[i]->AddMultMF(x, y);
370  }
371 }
372 
373 void SumIntegrator::AssembleDiagonalMF(Vector &diag)
374 {
375  for (int i = 0; i < integrators.Size(); i++)
376  {
377  integrators[i]->AssembleDiagonalMF(diag);
378  }
379 }
380 
381 void SumIntegrator::AssembleEA(const FiniteElementSpace &fes, Vector &emat,
382  const bool add)
383 {
384  for (int i = 0; i < integrators.Size(); i++)
385  {
386  integrators[i]->AssembleEA(fes, emat, add);
387  }
388 }
389 
390 void SumIntegrator::AssembleEAInteriorFaces(const FiniteElementSpace &fes,
391  Vector &ea_data_int,
392  Vector &ea_data_ext,
393  const bool add)
394 {
395  for (int i = 0; i < integrators.Size(); i++)
396  {
397  integrators[i]->AssembleEAInteriorFaces(fes,ea_data_int,ea_data_ext,add);
398  }
399 }
400 
401 void SumIntegrator::AssembleEABoundaryFaces(const FiniteElementSpace &fes,
402  Vector &ea_data_bdr,
403  const bool add)
404 {
405  for (int i = 0; i < integrators.Size(); i++)
406  {
407  integrators[i]->AssembleEABoundaryFaces(fes, ea_data_bdr, add);
408  }
409 }
410 
411 SumIntegrator::~SumIntegrator()
412 {
413  if (own_integrators)
414  {
415  for (int i = 0; i < integrators.Size(); i++)
416  {
417  delete integrators[i];
418  }
419  }
420 }
421 
422 void MixedScalarIntegrator::AssembleElementMatrix2(
423  const FiniteElement &trial_fe, const FiniteElement &test_fe,
425 {
426  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
427  this->FiniteElementTypeFailureMessage());
428 
429  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
430  bool same_shapes = same_calc_shape && (&trial_fe == &test_fe);
431 
432 #ifdef MFEM_THREAD_SAFE
433  Vector test_shape(test_nd);
434  Vector trial_shape;
435 #else
436  test_shape.SetSize(test_nd);
437 #endif
438  if (same_shapes)
439  {
440  trial_shape.NewDataAndSize(test_shape.GetData(), trial_nd);
441  }
442  else
443  {
444  trial_shape.SetSize(trial_nd);
445  }
446 
447  elmat.SetSize(test_nd, trial_nd);
448 
449  const IntegrationRule *ir = IntRule;
450  if (ir == NULL)
451  {
452  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
453  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
454  }
455 
456  elmat = 0.0;
457  for (i = 0; i < ir->GetNPoints(); i++)
458  {
459  const IntegrationPoint &ip = ir->IntPoint(i);
460  Trans.SetIntPoint(&ip);
461 
462  this->CalcTestShape(test_fe, Trans, test_shape);
463  this->CalcTrialShape(trial_fe, Trans, trial_shape);
464 
465  double w = Trans.Weight() * ip.weight;
466 
467  if (Q)
468  {
469  w *= Q->Eval(Trans, ip);
470  }
471  AddMult_a_VWt(w, test_shape, trial_shape, elmat);
472  }
473 #ifndef MFEM_THREAD_SAFE
474  if (same_shapes)
475  {
476  trial_shape.SetDataAndSize(NULL, 0);
477  }
478 #endif
479 }
480 
481 void MixedVectorIntegrator::AssembleElementMatrix2(
482  const FiniteElement &trial_fe, const FiniteElement &test_fe,
484 {
485  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
486  this->FiniteElementTypeFailureMessage());
487 
488  space_dim = Trans.GetSpaceDim();
489  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
490  int test_vdim = GetTestVDim(test_fe);
491  int trial_vdim = GetTrialVDim(trial_fe);
492  bool same_shapes = same_calc_shape && (&trial_fe == &test_fe);
493 
494  if (MQ)
495  {
496  MFEM_VERIFY(MQ->GetHeight() == test_vdim,
497  "Dimension mismatch in height of matrix coefficient.");
498  MFEM_VERIFY(MQ->GetWidth() == trial_vdim,
499  "Dimension mismatch in width of matrix coefficient.");
500  }
501  if (DQ)
502  {
503  MFEM_VERIFY(trial_vdim == test_vdim,
504  "Diagonal matrix coefficient requires matching "
505  "test and trial vector dimensions.");
506  MFEM_VERIFY(DQ->GetVDim() == trial_vdim,
507  "Dimension mismatch in diagonal matrix coefficient.");
508  }
509  if (VQ)
510  {
511  MFEM_VERIFY(VQ->GetVDim() == 3, "Vector coefficient must have "
512  "dimension equal to three.");
513  }
514 
515 #ifdef MFEM_THREAD_SAFE
516  Vector V(VQ ? VQ->GetVDim() : 0);
517  Vector D(DQ ? DQ->GetVDim() : 0);
518  DenseMatrix M(MQ ? MQ->GetHeight() : 0, MQ ? MQ->GetWidth() : 0);
519  DenseMatrix test_shape(test_nd, test_vdim);
520  DenseMatrix trial_shape;
521  DenseMatrix shape_tmp(test_nd, trial_vdim);
522 #else
523  V.SetSize(VQ ? VQ->GetVDim() : 0);
524  D.SetSize(DQ ? DQ->GetVDim() : 0);
525  M.SetSize(MQ ? MQ->GetHeight() : 0, MQ ? MQ->GetWidth() : 0);
526  test_shape.SetSize(test_nd, test_vdim);
527  shape_tmp.SetSize(test_nd, trial_vdim);
528 #endif
529  if (same_shapes)
530  {
531  trial_shape.Reset(test_shape.Data(), trial_nd, trial_vdim);
532  }
533  else
534  {
535  trial_shape.SetSize(trial_nd, trial_vdim);
536  }
537 
538  elmat.SetSize(test_nd, trial_nd);
539 
540  const IntegrationRule *ir = IntRule;
541  if (ir == NULL)
542  {
543  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
544  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
545  }
546 
547  elmat = 0.0;
548  for (i = 0; i < ir->GetNPoints(); i++)
549  {
550  const IntegrationPoint &ip = ir->IntPoint(i);
551  Trans.SetIntPoint(&ip);
552 
553  this->CalcTestShape(test_fe, Trans, test_shape);
554  if (!same_shapes)
555  {
556  this->CalcTrialShape(trial_fe, Trans, trial_shape);
557  }
558 
559  double w = Trans.Weight() * ip.weight;
560 
561  if (MQ)
562  {
563  MQ->Eval(M, Trans, ip);
564  M *= w;
565  Mult(test_shape, M, shape_tmp);
566  AddMultABt(shape_tmp, trial_shape, elmat);
567  }
568  else if (DQ)
569  {
570  DQ->Eval(D, Trans, ip);
571  D *= w;
572  AddMultADBt(test_shape, D, trial_shape, elmat);
573  }
574  else if (VQ)
575  {
576  VQ->Eval(V, Trans, ip);
577  V *= w;
578 
579  for (int j=0; j<test_nd; j++)
580  {
581  // Compute shape_tmp = test_shape x V
582  // V will always be of length 3
583  // shape_dim and test_shape could have reduced dimension
584  // i.e. 1D or 2D
585  if (test_vdim == 3 && trial_vdim == 3)
586  {
587  shape_tmp(j,0) = test_shape(j,1) * V(2) -
588  test_shape(j,2) * V(1);
589  shape_tmp(j,1) = test_shape(j,2) * V(0) -
590  test_shape(j,0) * V(2);
591  shape_tmp(j,2) = test_shape(j,0) * V(1) -
592  test_shape(j,1) * V(0);
593  }
594  else if (test_vdim == 3 && trial_vdim == 2)
595  {
596  shape_tmp(j,0) = test_shape(j,1) * V(2) -
597  test_shape(j,2) * V(1);
598  shape_tmp(j,1) = test_shape(j,2) * V(0) -
599  test_shape(j,0) * V(2);
600  }
601  else if (test_vdim == 3 && trial_vdim == 1)
602  {
603  shape_tmp(j,0) = test_shape(j,1) * V(2) -
604  test_shape(j,2) * V(1);
605  }
606  else if (test_vdim == 2 && trial_vdim == 3)
607  {
608  shape_tmp(j,0) = test_shape(j,1) * V(2);
609  shape_tmp(j,1) = -test_shape(j,0) * V(2);
610  shape_tmp(j,2) = test_shape(j,0) * V(1) -
611  test_shape(j,1) * V(0);
612  }
613  else if (test_vdim == 2 && trial_vdim == 2)
614  {
615  shape_tmp(j,0) = test_shape(j,1) * V(2);
616  shape_tmp(j,1) = -test_shape(j,0) * V(2);
617  }
618  else if (test_vdim == 1 && trial_vdim == 3)
619  {
620  shape_tmp(j,0) = 0.0;
621  shape_tmp(j,1) = -test_shape(j,0) * V(2);
622  shape_tmp(j,2) = test_shape(j,0) * V(1);
623  }
624  else if (test_vdim == 1 && trial_vdim == 1)
625  {
626  shape_tmp(j,0) = 0.0;
627  }
628  }
629  AddMultABt(shape_tmp, trial_shape, elmat);
630  }
631  else
632  {
633  if (Q)
634  {
635  w *= Q -> Eval (Trans, ip);
636  }
637  if (same_shapes)
638  {
639  AddMult_a_AAt (w, test_shape, elmat);
640  }
641  else
642  {
643  AddMult_a_ABt (w, test_shape, trial_shape, elmat);
644  }
645  }
646  }
647 #ifndef MFEM_THREAD_SAFE
648  if (same_shapes)
649  {
650  trial_shape.ClearExternalData();
651  }
652 #endif
653 }
654 
655 void MixedScalarVectorIntegrator::AssembleElementMatrix2(
656  const FiniteElement &trial_fe, const FiniteElement &test_fe,
658 {
659  MFEM_ASSERT(this->VerifyFiniteElementTypes(trial_fe, test_fe),
660  this->FiniteElementTypeFailureMessage());
661 
662  MFEM_VERIFY(VQ, "MixedScalarVectorIntegrator: "
663  "VectorCoefficient must be set");
664 
665  const FiniteElement * vec_fe = transpose?&trial_fe:&test_fe;
666  const FiniteElement * sca_fe = transpose?&test_fe:&trial_fe;
667 
668  space_dim = Trans.GetSpaceDim();
669  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
670  int sca_nd = sca_fe->GetDof();
671  int vec_nd = vec_fe->GetDof();
672  int vdim = GetVDim(*vec_fe);
673  double vtmp;
674 
675  MFEM_VERIFY(VQ->GetVDim() == vdim, "MixedScalarVectorIntegrator: "
676  "Dimensions of VectorCoefficient and Vector-valued basis "
677  "functions must match");
678 
679 #ifdef MFEM_THREAD_SAFE
680  Vector V(vdim);
681  DenseMatrix vshape(vec_nd, vdim);
682  Vector shape(sca_nd);
683  Vector vshape_tmp(vec_nd);
684 #else
685  V.SetSize(vdim);
686  vshape.SetSize(vec_nd, vdim);
687  shape.SetSize(sca_nd);
688  vshape_tmp.SetSize(vec_nd);
689 #endif
690 
691  Vector V_test(transpose?shape.GetData():vshape_tmp.GetData(),test_nd);
692  Vector W_trial(transpose?vshape_tmp.GetData():shape.GetData(),trial_nd);
693 
694  elmat.SetSize(test_nd, trial_nd);
695 
696  const IntegrationRule *ir = IntRule;
697  if (ir == NULL)
698  {
699  int ir_order = this->GetIntegrationOrder(trial_fe, test_fe, Trans);
700  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
701  }
702 
703  elmat = 0.0;
704  for (i = 0; i < ir->GetNPoints(); i++)
705  {
706  const IntegrationPoint &ip = ir->IntPoint(i);
707  Trans.SetIntPoint(&ip);
708 
709  this->CalcShape(*sca_fe, Trans, shape);
710  this->CalcVShape(*vec_fe, Trans, vshape);
711 
712  double w = Trans.Weight() * ip.weight;
713 
714  VQ->Eval(V, Trans, ip);
715  V *= w;
716 
717  if ( vdim == 2 && cross_2d )
718  {
719  vtmp = V[0];
720  V[0] = -V[1];
721  V[1] = vtmp;
722  }
723 
724  vshape.Mult(V,vshape_tmp);
725  AddMultVWt(V_test, W_trial, elmat);
726  }
727 }
728 
729 
730 void GradientIntegrator::AssembleElementMatrix2(
731  const FiniteElement &trial_fe, const FiniteElement &test_fe,
733 {
734  dim = test_fe.GetDim();
735  int trial_dof = trial_fe.GetDof();
736  int test_dof = test_fe.GetDof();
737  double c;
738  Vector d_col;
739 
740  dshape.SetSize(trial_dof, dim);
741  gshape.SetSize(trial_dof, dim);
742  Jadj.SetSize(dim);
743  shape.SetSize(test_dof);
744  elmat.SetSize(dim * test_dof, trial_dof);
745 
746  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(trial_fe, test_fe,
747  Trans);
748 
749  elmat = 0.0;
750  elmat_comp.SetSize(test_dof, trial_dof);
751 
752  for (int i = 0; i < ir->GetNPoints(); i++)
753  {
754  const IntegrationPoint &ip = ir->IntPoint(i);
755 
756  trial_fe.CalcDShape(ip, dshape);
757  test_fe.CalcShape(ip, shape);
758 
759  Trans.SetIntPoint(&ip);
760  CalcAdjugate(Trans.Jacobian(), Jadj);
761 
762  Mult(dshape, Jadj, gshape);
763 
764  c = ip.weight;
765  if (Q)
766  {
767  c *= Q->Eval(Trans, ip);
768  }
769  shape *= c;
770 
771  for (int d = 0; d < dim; ++d)
772  {
773  gshape.GetColumnReference(d, d_col);
774  MultVWt(shape, d_col, elmat_comp);
775  for (int jj = 0; jj < trial_dof; ++jj)
776  {
777  for (int ii = 0; ii < test_dof; ++ii)
778  {
779  elmat(d * test_dof + ii, jj) += elmat_comp(ii, jj);
780  }
781  }
782  }
783  }
784 }
785 
786 const IntegrationRule &GradientIntegrator::GetRule(const FiniteElement
787  &trial_fe,
788  const FiniteElement &test_fe,
790 {
791  int order = Trans.OrderGrad(&trial_fe) + test_fe.GetOrder() + Trans.OrderJ();
792  return IntRules.Get(trial_fe.GetGeomType(), order);
793 }
794 
795 
796 void DiffusionIntegrator::AssembleElementMatrix
798  DenseMatrix &elmat )
799 {
800  int nd = el.GetDof();
801  dim = el.GetDim();
802  int spaceDim = Trans.GetSpaceDim();
803  bool square = (dim == spaceDim);
804  double w;
805 
806  if (VQ)
807  {
808  MFEM_VERIFY(VQ->GetVDim() == spaceDim,
809  "Unexpected dimension for VectorCoefficient");
810  }
811  if (MQ)
812  {
813  MFEM_VERIFY(MQ->GetWidth() == spaceDim,
814  "Unexpected width for MatrixCoefficient");
815  MFEM_VERIFY(MQ->GetHeight() == spaceDim,
816  "Unexpected height for MatrixCoefficient");
817  }
818 
819 #ifdef MFEM_THREAD_SAFE
820  DenseMatrix dshape(nd, dim), dshapedxt(nd, spaceDim);
821  DenseMatrix dshapedxt_m(nd, MQ ? spaceDim : 0);
822  DenseMatrix M(MQ ? spaceDim : 0);
823  Vector D(VQ ? VQ->GetVDim() : 0);
824 #else
825  dshape.SetSize(nd, dim);
826  dshapedxt.SetSize(nd, spaceDim);
827  dshapedxt_m.SetSize(nd, MQ ? spaceDim : 0);
828  M.SetSize(MQ ? spaceDim : 0);
829  D.SetSize(VQ ? VQ->GetVDim() : 0);
830 #endif
831  elmat.SetSize(nd);
832 
833  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(el, el);
834 
835  elmat = 0.0;
836  for (int i = 0; i < ir->GetNPoints(); i++)
837  {
838  const IntegrationPoint &ip = ir->IntPoint(i);
839  el.CalcDShape(ip, dshape);
840 
841  Trans.SetIntPoint(&ip);
842  w = Trans.Weight();
843  w = ip.weight / (square ? w : w*w*w);
844  // AdjugateJacobian = / adj(J), if J is square
845  // \ adj(J^t.J).J^t, otherwise
846  Mult(dshape, Trans.AdjugateJacobian(), dshapedxt);
847  if (MQ)
848  {
849  MQ->Eval(M, Trans, ip);
850  M *= w;
851  Mult(dshapedxt, M, dshapedxt_m);
852  AddMultABt(dshapedxt_m, dshapedxt, elmat);
853  }
854  else if (VQ)
855  {
856  VQ->Eval(D, Trans, ip);
857  D *= w;
858  AddMultADAt(dshapedxt, D, elmat);
859  }
860  else
861  {
862  if (Q)
863  {
864  w *= Q->Eval(Trans, ip);
865  }
866  AddMult_a_AAt(w, dshapedxt, elmat);
867  }
868  }
869 }
870 
871 void DiffusionIntegrator::AssembleElementMatrix2(
872  const FiniteElement &trial_fe, const FiniteElement &test_fe,
873  ElementTransformation &Trans, DenseMatrix &elmat)
874 {
875  int tr_nd = trial_fe.GetDof();
876  int te_nd = test_fe.GetDof();
877  dim = trial_fe.GetDim();
878  int spaceDim = Trans.GetSpaceDim();
879  bool square = (dim == spaceDim);
880  double w;
881 
882  if (VQ)
883  {
884  MFEM_VERIFY(VQ->GetVDim() == spaceDim,
885  "Unexpected dimension for VectorCoefficient");
886  }
887  if (MQ)
888  {
889  MFEM_VERIFY(MQ->GetWidth() == spaceDim,
890  "Unexpected width for MatrixCoefficient");
891  MFEM_VERIFY(MQ->GetHeight() == spaceDim,
892  "Unexpected height for MatrixCoefficient");
893  }
894 
895 #ifdef MFEM_THREAD_SAFE
896  DenseMatrix dshape(tr_nd, dim), dshapedxt(tr_nd, spaceDim);
897  DenseMatrix te_dshape(te_nd, dim), te_dshapedxt(te_nd, spaceDim);
898  DenseMatrix invdfdx(dim, spaceDim);
899  DenseMatrix dshapedxt_m(te_nd, MQ ? spaceDim : 0);
900  DenseMatrix M(MQ ? spaceDim : 0);
901  Vector D(VQ ? VQ->GetVDim() : 0);
902 #else
903  dshape.SetSize(tr_nd, dim);
904  dshapedxt.SetSize(tr_nd, spaceDim);
905  te_dshape.SetSize(te_nd, dim);
906  te_dshapedxt.SetSize(te_nd, spaceDim);
907  invdfdx.SetSize(dim, spaceDim);
908  dshapedxt_m.SetSize(te_nd, MQ ? spaceDim : 0);
909  M.SetSize(MQ ? spaceDim : 0);
910  D.SetSize(VQ ? VQ->GetVDim() : 0);
911 #endif
912  elmat.SetSize(te_nd, tr_nd);
913 
914  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(trial_fe, test_fe);
915 
916  elmat = 0.0;
917  for (int i = 0; i < ir->GetNPoints(); i++)
918  {
919  const IntegrationPoint &ip = ir->IntPoint(i);
920  trial_fe.CalcDShape(ip, dshape);
921  test_fe.CalcDShape(ip, te_dshape);
922 
923  Trans.SetIntPoint(&ip);
924  CalcAdjugate(Trans.Jacobian(), invdfdx);
925  w = Trans.Weight();
926  w = ip.weight / (square ? w : w*w*w);
927  Mult(dshape, invdfdx, dshapedxt);
928  Mult(te_dshape, invdfdx, te_dshapedxt);
929  // invdfdx, dshape, and te_dshape no longer needed
930  if (MQ)
931  {
932  MQ->Eval(M, Trans, ip);
933  M *= w;
934  Mult(te_dshapedxt, M, dshapedxt_m);
935  AddMultABt(dshapedxt_m, dshapedxt, elmat);
936  }
937  else if (VQ)
938  {
939  VQ->Eval(D, Trans, ip);
940  D *= w;
941  AddMultADAt(dshapedxt, D, elmat);
942  }
943  else
944  {
945  if (Q)
946  {
947  w *= Q->Eval(Trans, ip);
948  }
949  dshapedxt *= w;
950  AddMultABt(te_dshapedxt, dshapedxt, elmat);
951  }
952  }
953 }
954 
955 void DiffusionIntegrator::AssembleElementVector(
956  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
957  Vector &elvect)
958 {
959  int nd = el.GetDof();
960  dim = el.GetDim();
961  int spaceDim = Tr.GetSpaceDim();
962  double w;
963 
964  if (VQ)
965  {
966  MFEM_VERIFY(VQ->GetVDim() == spaceDim,
967  "Unexpected dimension for VectorCoefficient");
968  }
969  if (MQ)
970  {
971  MFEM_VERIFY(MQ->GetWidth() == spaceDim,
972  "Unexpected width for MatrixCoefficient");
973  MFEM_VERIFY(MQ->GetHeight() == spaceDim,
974  "Unexpected height for MatrixCoefficient");
975  }
976 
977 #ifdef MFEM_THREAD_SAFE
978  DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim), M(MQ ? spaceDim : 0);
979  Vector D(VQ ? VQ->GetVDim() : 0);
980 #else
981  dshape.SetSize(nd,dim);
982  invdfdx.SetSize(dim, spaceDim);
983  M.SetSize(MQ ? spaceDim : 0);
984  D.SetSize(VQ ? VQ->GetVDim() : 0);
985 #endif
986  vec.SetSize(dim);
987  vecdxt.SetSize((VQ || MQ) ? spaceDim : 0);
988  pointflux.SetSize(spaceDim);
989 
990  elvect.SetSize(nd);
991 
992  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(el, el);
993 
994  elvect = 0.0;
995  for (int i = 0; i < ir->GetNPoints(); i++)
996  {
997  const IntegrationPoint &ip = ir->IntPoint(i);
998  el.CalcDShape(ip, dshape);
999 
1000  Tr.SetIntPoint(&ip);
1001  CalcAdjugate(Tr.Jacobian(), invdfdx); // invdfdx = adj(J)
1002  w = ip.weight / Tr.Weight();
1003 
1004  if (!MQ && !VQ)
1005  {
1006  dshape.MultTranspose(elfun, vec);
1007  invdfdx.MultTranspose(vec, pointflux);
1008  if (Q)
1009  {
1010  w *= Q->Eval(Tr, ip);
1011  }
1012  }
1013  else
1014  {
1015  dshape.MultTranspose(elfun, vec);
1016  invdfdx.MultTranspose(vec, vecdxt);
1017  if (MQ)
1018  {
1019  MQ->Eval(M, Tr, ip);
1020  M.Mult(vecdxt, pointflux);
1021  }
1022  else
1023  {
1024  VQ->Eval(D, Tr, ip);
1025  for (int j=0; j<spaceDim; ++j)
1026  {
1027  pointflux[j] = D[j] * vecdxt[j];
1028  }
1029  }
1030  }
1031  pointflux *= w;
1032  invdfdx.Mult(pointflux, vec);
1033  dshape.AddMult(vec, elvect);
1034  }
1035 }
1036 
1037 void DiffusionIntegrator::ComputeElementFlux
1039  Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef )
1040 {
1041  int nd, spaceDim, fnd;
1042 
1043  nd = el.GetDof();
1044  dim = el.GetDim();
1045  spaceDim = Trans.GetSpaceDim();
1046 
1047  if (VQ)
1048  {
1049  MFEM_VERIFY(VQ->GetVDim() == spaceDim,
1050  "Unexpected dimension for VectorCoefficient");
1051  }
1052  if (MQ)
1053  {
1054  MFEM_VERIFY(MQ->GetWidth() == spaceDim,
1055  "Unexpected width for MatrixCoefficient");
1056  MFEM_VERIFY(MQ->GetHeight() == spaceDim,
1057  "Unexpected height for MatrixCoefficient");
1058  }
1059 
1060  MFEM_VERIFY(!SMQ, "SymmetricMatrixCoefficient not supported here");
1061 
1062 #ifdef MFEM_THREAD_SAFE
1063  DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim);
1064  DenseMatrix M(MQ ? spaceDim : 0);
1065  Vector D(VQ ? VQ->GetVDim() : 0);
1066 #else
1067  dshape.SetSize(nd,dim);
1068  invdfdx.SetSize(dim, spaceDim);
1069  M.SetSize(MQ ? spaceDim : 0);
1070  D.SetSize(VQ ? VQ->GetVDim() : 0);
1071 #endif
1072  vec.SetSize(dim);
1073  vecdxt.SetSize(spaceDim);
1074  pointflux.SetSize(MQ || VQ ? spaceDim : 0);
1075 
1076  const IntegrationRule &ir = fluxelem.GetNodes();
1077  fnd = ir.GetNPoints();
1078  flux.SetSize( fnd * spaceDim );
1079 
1080  for (int i = 0; i < fnd; i++)
1081  {
1082  const IntegrationPoint &ip = ir.IntPoint(i);
1083  el.CalcDShape(ip, dshape);
1084  dshape.MultTranspose(u, vec);
1085 
1086  Trans.SetIntPoint (&ip);
1087  CalcInverse(Trans.Jacobian(), invdfdx);
1088  invdfdx.MultTranspose(vec, vecdxt);
1089 
1090  if (with_coef)
1091  {
1092  if (!MQ && !VQ)
1093  {
1094  if (Q)
1095  {
1096  vecdxt *= Q->Eval(Trans,ip);
1097  }
1098  for (int j = 0; j < spaceDim; j++)
1099  {
1100  flux(fnd*j+i) = vecdxt(j);
1101  }
1102  }
1103  else
1104  {
1105  if (MQ)
1106  {
1107  MQ->Eval(M, Trans, ip);
1108  M.Mult(vecdxt, pointflux);
1109  }
1110  else
1111  {
1112  VQ->Eval(D, Trans, ip);
1113  for (int j=0; j<spaceDim; ++j)
1114  {
1115  pointflux[j] = D[j] * vecdxt[j];
1116  }
1117  }
1118  for (int j = 0; j < spaceDim; j++)
1119  {
1120  flux(fnd*j+i) = pointflux(j);
1121  }
1122  }
1123  }
1124  else
1125  {
1126  for (int j = 0; j < spaceDim; j++)
1127  {
1128  flux(fnd*j+i) = vecdxt(j);
1129  }
1130  }
1131  }
1132 }
1133 
1134 double DiffusionIntegrator::ComputeFluxEnergy
1135 ( const FiniteElement &fluxelem, ElementTransformation &Trans,
1136  Vector &flux, Vector* d_energy)
1137 {
1138  int nd = fluxelem.GetDof();
1139  dim = fluxelem.GetDim();
1140  int spaceDim = Trans.GetSpaceDim();
1141 
1142 #ifdef MFEM_THREAD_SAFE
1143  DenseMatrix M;
1144  Vector D(VQ ? VQ->GetVDim() : 0);
1145 #else
1146  D.SetSize(VQ ? VQ->GetVDim() : 0);
1147 #endif
1148 
1149  MFEM_VERIFY(!SMQ, "SymmetricMatrixCoefficient not supported here");
1150 
1151  shape.SetSize(nd);
1152  pointflux.SetSize(spaceDim);
1153  if (d_energy) { vec.SetSize(spaceDim); }
1154  if (MQ) { M.SetSize(spaceDim); }
1155 
1156  int order = 2 * fluxelem.GetOrder(); // <--
1157  const IntegrationRule *ir = &IntRules.Get(fluxelem.GetGeomType(), order);
1158 
1159  double energy = 0.0;
1160  if (d_energy) { *d_energy = 0.0; }
1161 
1162  for (int i = 0; i < ir->GetNPoints(); i++)
1163  {
1164  const IntegrationPoint &ip = ir->IntPoint(i);
1165  fluxelem.CalcShape(ip, shape);
1166 
1167  pointflux = 0.0;
1168  for (int k = 0; k < spaceDim; k++)
1169  {
1170  for (int j = 0; j < nd; j++)
1171  {
1172  pointflux(k) += flux(k*nd+j)*shape(j);
1173  }
1174  }
1175 
1176  Trans.SetIntPoint(&ip);
1177  double w = Trans.Weight() * ip.weight;
1178 
1179  if (MQ)
1180  {
1181  MQ->Eval(M, Trans, ip);
1182  energy += w * M.InnerProduct(pointflux, pointflux);
1183  }
1184  else if (VQ)
1185  {
1186  VQ->Eval(D, Trans, ip);
1187  D *= pointflux;
1188  energy += w * (D * pointflux);
1189  }
1190  else
1191  {
1192  double e = (pointflux * pointflux);
1193  if (Q) { e *= Q->Eval(Trans, ip); }
1194  energy += w * e;
1195  }
1196 
1197  if (d_energy)
1198  {
1199  // transform pointflux to the ref. domain and integrate the components
1200  Trans.Jacobian().MultTranspose(pointflux, vec);
1201  for (int k = 0; k < dim; k++)
1202  {
1203  (*d_energy)[k] += w * vec[k] * vec[k];
1204  }
1205  // TODO: Q, VQ, MQ
1206  }
1207  }
1208 
1209  return energy;
1210 }
1211 
1212 const IntegrationRule &DiffusionIntegrator::GetRule(
1213  const FiniteElement &trial_fe, const FiniteElement &test_fe)
1214 {
1215  int order;
1216  if (trial_fe.Space() == FunctionSpace::Pk)
1217  {
1218  order = trial_fe.GetOrder() + test_fe.GetOrder() - 2;
1219  }
1220  else
1221  {
1222  // order = 2*el.GetOrder() - 2; // <-- this seems to work fine too
1223  order = trial_fe.GetOrder() + test_fe.GetOrder() + trial_fe.GetDim() - 1;
1224  }
1225 
1226  if (trial_fe.Space() == FunctionSpace::rQk)
1227  {
1228  return RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1229  }
1230  return IntRules.Get(trial_fe.GetGeomType(), order);
1231 }
1232 
1233 
1234 void MassIntegrator::AssembleElementMatrix
1236  DenseMatrix &elmat )
1237 {
1238  int nd = el.GetDof();
1239  // int dim = el.GetDim();
1240  double w;
1241 
1242 #ifdef MFEM_THREAD_SAFE
1243  Vector shape;
1244 #endif
1245  elmat.SetSize(nd);
1246  shape.SetSize(nd);
1247 
1248  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(el, el, Trans);
1249 
1250  elmat = 0.0;
1251  for (int i = 0; i < ir->GetNPoints(); i++)
1252  {
1253  const IntegrationPoint &ip = ir->IntPoint(i);
1254  Trans.SetIntPoint (&ip);
1255 
1256  el.CalcPhysShape(Trans, shape);
1257 
1258  w = Trans.Weight() * ip.weight;
1259  if (Q)
1260  {
1261  w *= Q -> Eval(Trans, ip);
1262  }
1263 
1264  AddMult_a_VVt(w, shape, elmat);
1265  }
1266 }
1267 
1268 void MassIntegrator::AssembleElementMatrix2(
1269  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1270  ElementTransformation &Trans, DenseMatrix &elmat)
1271 {
1272  int tr_nd = trial_fe.GetDof();
1273  int te_nd = test_fe.GetDof();
1274  double w;
1275 
1276 #ifdef MFEM_THREAD_SAFE
1277  Vector shape, te_shape;
1278 #endif
1279  elmat.SetSize(te_nd, tr_nd);
1280  shape.SetSize(tr_nd);
1281  te_shape.SetSize(te_nd);
1282 
1283  const IntegrationRule *ir = IntRule ? IntRule :
1284  &GetRule(trial_fe, test_fe, Trans);
1285 
1286  elmat = 0.0;
1287  for (int i = 0; i < ir->GetNPoints(); i++)
1288  {
1289  const IntegrationPoint &ip = ir->IntPoint(i);
1290  trial_fe.CalcShape(ip, shape);
1291  test_fe.CalcShape(ip, te_shape);
1292 
1293  Trans.SetIntPoint (&ip);
1294  w = Trans.Weight() * ip.weight;
1295  if (Q)
1296  {
1297  w *= Q -> Eval(Trans, ip);
1298  }
1299 
1300  te_shape *= w;
1301  AddMultVWt(te_shape, shape, elmat);
1302  }
1303 }
1304 
1305 const IntegrationRule &MassIntegrator::GetRule(const FiniteElement &trial_fe,
1306  const FiniteElement &test_fe,
1307  ElementTransformation &Trans)
1308 {
1309  // int order = trial_fe.GetOrder() + test_fe.GetOrder();
1310  const int order = trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW();
1311 
1312  if (trial_fe.Space() == FunctionSpace::rQk)
1313  {
1314  return RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1315  }
1316  return IntRules.Get(trial_fe.GetGeomType(), order);
1317 }
1318 
1319 
1320 void BoundaryMassIntegrator::AssembleFaceMatrix(
1321  const FiniteElement &el1, const FiniteElement &el2,
1322  FaceElementTransformations &Trans, DenseMatrix &elmat)
1323 {
1324  MFEM_ASSERT(Trans.Elem2No < 0,
1325  "support for interior faces is not implemented");
1326 
1327  int nd1 = el1.GetDof();
1328  double w;
1329 
1330 #ifdef MFEM_THREAD_SAFE
1331  Vector shape;
1332 #endif
1333  elmat.SetSize(nd1);
1334  shape.SetSize(nd1);
1335 
1336  const IntegrationRule *ir = IntRule;
1337  if (ir == NULL)
1338  {
1339  int order = 2 * el1.GetOrder();
1340 
1341  ir = &IntRules.Get(Trans.GetGeometryType(), order);
1342  }
1343 
1344  elmat = 0.0;
1345  for (int i = 0; i < ir->GetNPoints(); i++)
1346  {
1347  const IntegrationPoint &ip = ir->IntPoint(i);
1348 
1349  // Set the integration point in the face and the neighboring element
1350  Trans.SetAllIntPoints(&ip);
1351 
1352  // Access the neighboring element's integration point
1353  const IntegrationPoint &eip = Trans.GetElement1IntPoint();
1354  el1.CalcShape(eip, shape);
1355 
1356  w = Trans.Weight() * ip.weight;
1357  if (Q)
1358  {
1359  w *= Q -> Eval(Trans, ip);
1360  }
1361 
1362  AddMult_a_VVt(w, shape, elmat);
1363  }
1364 }
1365 
1366 void ConvectionIntegrator::AssembleElementMatrix(
1367  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1368 {
1369  int nd = el.GetDof();
1370  dim = el.GetDim();
1371 
1372 #ifdef MFEM_THREAD_SAFE
1373  DenseMatrix dshape, adjJ, Q_ir;
1374  Vector shape, vec2, BdFidxT;
1375 #endif
1376  elmat.SetSize(nd);
1377  dshape.SetSize(nd,dim);
1378  adjJ.SetSize(dim);
1379  shape.SetSize(nd);
1380  vec2.SetSize(dim);
1381  BdFidxT.SetSize(nd);
1382 
1383  Vector vec1;
1384 
1385  const IntegrationRule *ir = IntRule;
1386  if (ir == NULL)
1387  {
1388  int order = Trans.OrderGrad(&el) + Trans.Order() + el.GetOrder();
1389  ir = &IntRules.Get(el.GetGeomType(), order);
1390  }
1391 
1392  Q->Eval(Q_ir, Trans, *ir);
1393 
1394  elmat = 0.0;
1395  for (int i = 0; i < ir->GetNPoints(); i++)
1396  {
1397  const IntegrationPoint &ip = ir->IntPoint(i);
1398  el.CalcDShape(ip, dshape);
1399  el.CalcShape(ip, shape);
1400 
1401  Trans.SetIntPoint(&ip);
1402  CalcAdjugate(Trans.Jacobian(), adjJ);
1403  Q_ir.GetColumnReference(i, vec1);
1404  vec1 *= alpha * ip.weight;
1405 
1406  adjJ.Mult(vec1, vec2);
1407  dshape.Mult(vec2, BdFidxT);
1408 
1409  AddMultVWt(shape, BdFidxT, elmat);
1410  }
1411 }
1412 
1413 
1414 void GroupConvectionIntegrator::AssembleElementMatrix(
1415  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1416 {
1417  int nd = el.GetDof();
1418  int dim = el.GetDim();
1419 
1420  elmat.SetSize(nd);
1421  dshape.SetSize(nd,dim);
1422  adjJ.SetSize(dim);
1423  shape.SetSize(nd);
1424  grad.SetSize(nd,dim);
1425 
1426  const IntegrationRule *ir = IntRule;
1427  if (ir == NULL)
1428  {
1429  int order = Trans.OrderGrad(&el) + el.GetOrder();
1430  ir = &IntRules.Get(el.GetGeomType(), order);
1431  }
1432 
1433  Q->Eval(Q_nodal, Trans, el.GetNodes()); // sets the size of Q_nodal
1434 
1435  elmat = 0.0;
1436  for (int i = 0; i < ir->GetNPoints(); i++)
1437  {
1438  const IntegrationPoint &ip = ir->IntPoint(i);
1439  el.CalcDShape(ip, dshape);
1440  el.CalcShape(ip, shape);
1441 
1442  Trans.SetIntPoint(&ip);
1443  CalcAdjugate(Trans.Jacobian(), adjJ);
1444 
1445  Mult(dshape, adjJ, grad);
1446 
1447  double w = alpha * ip.weight;
1448 
1449  // elmat(k,l) += \sum_s w*shape(k)*Q_nodal(s,k)*grad(l,s)
1450  for (int k = 0; k < nd; k++)
1451  {
1452  double wsk = w*shape(k);
1453  for (int l = 0; l < nd; l++)
1454  {
1455  double a = 0.0;
1456  for (int s = 0; s < dim; s++)
1457  {
1458  a += Q_nodal(s,k)*grad(l,s);
1459  }
1460  elmat(k,l) += wsk*a;
1461  }
1462  }
1463  }
1464 }
1465 
1466 const IntegrationRule &ConvectionIntegrator::GetRule(
1467  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1468  ElementTransformation &Trans)
1469 {
1470  int order = Trans.OrderGrad(&trial_fe) + Trans.Order() + test_fe.GetOrder();
1471 
1472  return IntRules.Get(trial_fe.GetGeomType(), order);
1473 }
1474 
1475 const IntegrationRule &ConvectionIntegrator::GetRule(
1476  const FiniteElement &el, ElementTransformation &Trans)
1477 {
1478  return GetRule(el,el,Trans);
1479 }
1480 
1481 void VectorMassIntegrator::AssembleElementMatrix
1483  DenseMatrix &elmat )
1484 {
1485  int nd = el.GetDof();
1486  int spaceDim = Trans.GetSpaceDim();
1487 
1488  double norm;
1489 
1490  // If vdim is not set, set it to the space dimension
1491  vdim = (vdim == -1) ? spaceDim : vdim;
1492 
1493  elmat.SetSize(nd*vdim);
1494  shape.SetSize(nd);
1495  partelmat.SetSize(nd);
1496  if (VQ)
1497  {
1498  vec.SetSize(vdim);
1499  }
1500  else if (MQ)
1501  {
1502  mcoeff.SetSize(vdim);
1503  }
1504 
1505  const IntegrationRule *ir = IntRule;
1506  if (ir == NULL)
1507  {
1508  int order = 2 * el.GetOrder() + Trans.OrderW() + Q_order;
1509 
1510  if (el.Space() == FunctionSpace::rQk)
1511  {
1512  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
1513  }
1514  else
1515  {
1516  ir = &IntRules.Get(el.GetGeomType(), order);
1517  }
1518  }
1519 
1520  elmat = 0.0;
1521  for (int s = 0; s < ir->GetNPoints(); s++)
1522  {
1523  const IntegrationPoint &ip = ir->IntPoint(s);
1524  el.CalcShape(ip, shape);
1525 
1526  Trans.SetIntPoint (&ip);
1527  norm = ip.weight * Trans.Weight();
1528 
1529  MultVVt(shape, partelmat);
1530 
1531  if (VQ)
1532  {
1533  VQ->Eval(vec, Trans, ip);
1534  for (int k = 0; k < vdim; k++)
1535  {
1536  elmat.AddMatrix(norm*vec(k), partelmat, nd*k, nd*k);
1537  }
1538  }
1539  else if (MQ)
1540  {
1541  MQ->Eval(mcoeff, Trans, ip);
1542  for (int i = 0; i < vdim; i++)
1543  for (int j = 0; j < vdim; j++)
1544  {
1545  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, nd*i, nd*j);
1546  }
1547  }
1548  else
1549  {
1550  if (Q)
1551  {
1552  norm *= Q->Eval(Trans, ip);
1553  }
1554  partelmat *= norm;
1555  for (int k = 0; k < vdim; k++)
1556  {
1557  elmat.AddMatrix(partelmat, nd*k, nd*k);
1558  }
1559  }
1560  }
1561 }
1562 
1563 void VectorMassIntegrator::AssembleElementMatrix2(
1564  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1565  ElementTransformation &Trans, DenseMatrix &elmat)
1566 {
1567  int tr_nd = trial_fe.GetDof();
1568  int te_nd = test_fe.GetDof();
1569 
1570  double norm;
1571 
1572  // If vdim is not set, set it to the space dimension
1573  vdim = (vdim == -1) ? Trans.GetSpaceDim() : vdim;
1574 
1575  elmat.SetSize(te_nd*vdim, tr_nd*vdim);
1576  shape.SetSize(tr_nd);
1577  te_shape.SetSize(te_nd);
1578  partelmat.SetSize(te_nd, tr_nd);
1579  if (VQ)
1580  {
1581  vec.SetSize(vdim);
1582  }
1583  else if (MQ)
1584  {
1585  mcoeff.SetSize(vdim);
1586  }
1587 
1588  const IntegrationRule *ir = IntRule;
1589  if (ir == NULL)
1590  {
1591  int order = (trial_fe.GetOrder() + test_fe.GetOrder() +
1592  Trans.OrderW() + Q_order);
1593 
1594  if (trial_fe.Space() == FunctionSpace::rQk)
1595  {
1596  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1597  }
1598  else
1599  {
1600  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1601  }
1602  }
1603 
1604  elmat = 0.0;
1605  for (int s = 0; s < ir->GetNPoints(); s++)
1606  {
1607  const IntegrationPoint &ip = ir->IntPoint(s);
1608  trial_fe.CalcShape(ip, shape);
1609  test_fe.CalcShape(ip, te_shape);
1610 
1611  Trans.SetIntPoint(&ip);
1612  norm = ip.weight * Trans.Weight();
1613 
1614  MultVWt(te_shape, shape, partelmat);
1615 
1616  if (VQ)
1617  {
1618  VQ->Eval(vec, Trans, ip);
1619  for (int k = 0; k < vdim; k++)
1620  {
1621  elmat.AddMatrix(norm*vec(k), partelmat, te_nd*k, tr_nd*k);
1622  }
1623  }
1624  else if (MQ)
1625  {
1626  MQ->Eval(mcoeff, Trans, ip);
1627  for (int i = 0; i < vdim; i++)
1628  for (int j = 0; j < vdim; j++)
1629  {
1630  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, te_nd*i, tr_nd*j);
1631  }
1632  }
1633  else
1634  {
1635  if (Q)
1636  {
1637  norm *= Q->Eval(Trans, ip);
1638  }
1639  partelmat *= norm;
1640  for (int k = 0; k < vdim; k++)
1641  {
1642  elmat.AddMatrix(partelmat, te_nd*k, tr_nd*k);
1643  }
1644  }
1645  }
1646 }
1647 
1648 void VectorFEDivergenceIntegrator::AssembleElementMatrix2(
1649  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1650  ElementTransformation &Trans, DenseMatrix &elmat)
1651 {
1652  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1653 
1654 #ifdef MFEM_THREAD_SAFE
1655  Vector divshape(trial_nd), shape(test_nd);
1656 #else
1657  divshape.SetSize(trial_nd);
1658  shape.SetSize(test_nd);
1659 #endif
1660 
1661  elmat.SetSize(test_nd, trial_nd);
1662 
1663  const IntegrationRule *ir = IntRule;
1664  if (ir == NULL)
1665  {
1666  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
1667  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1668  }
1669 
1670  elmat = 0.0;
1671  for (i = 0; i < ir->GetNPoints(); i++)
1672  {
1673  const IntegrationPoint &ip = ir->IntPoint(i);
1674  trial_fe.CalcDivShape(ip, divshape);
1675  test_fe.CalcShape(ip, shape);
1676  double w = ip.weight;
1677  if (Q)
1678  {
1679  Trans.SetIntPoint(&ip);
1680  w *= Q->Eval(Trans, ip);
1681  }
1682  shape *= w;
1683  AddMultVWt(shape, divshape, elmat);
1684  }
1685 }
1686 
1687 void VectorFEWeakDivergenceIntegrator::AssembleElementMatrix2(
1688  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1689  ElementTransformation &Trans, DenseMatrix &elmat)
1690 {
1691  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1692  int dim = trial_fe.GetDim();
1693 
1694  MFEM_ASSERT(test_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1695  test_fe.GetMapType() == mfem::FiniteElement::VALUE &&
1697  "Trial space must be H(Curl) and test space must be H_1");
1698 
1699 #ifdef MFEM_THREAD_SAFE
1700  DenseMatrix dshape(test_nd, dim);
1701  DenseMatrix dshapedxt(test_nd, dim);
1702  DenseMatrix vshape(trial_nd, dim);
1703  DenseMatrix invdfdx(dim);
1704 #else
1705  dshape.SetSize(test_nd, dim);
1706  dshapedxt.SetSize(test_nd, dim);
1707  vshape.SetSize(trial_nd, dim);
1708  invdfdx.SetSize(dim);
1709 #endif
1710 
1711  elmat.SetSize(test_nd, trial_nd);
1712 
1713  const IntegrationRule *ir = IntRule;
1714  if (ir == NULL)
1715  {
1716  // The integrand on the reference element is:
1717  // -( Q/det(J) ) u_hat^T adj(J) adj(J)^T grad_hat(v_hat).
1718  //
1719  // For Trans in (P_k)^d, v_hat in P_l, u_hat in ND_m, and dim=sdim=d>=1
1720  // - J_{ij} is in P_{k-1}, so adj(J)_{ij} is in P_{(d-1)*(k-1)}
1721  // - so adj(J)^T grad_hat(v_hat) is in (P_{(d-1)*(k-1)+(l-1)})^d
1722  // - u_hat is in (P_m)^d
1723  // - adj(J)^T u_hat is in (P_{(d-1)*(k-1)+m})^d
1724  // - and u_hat^T adj(J) adj(J)^T grad_hat(v_hat) is in P_n with
1725  // n = 2*(d-1)*(k-1)+(l-1)+m
1726  //
1727  // For Trans in (Q_k)^d, v_hat in Q_l, u_hat in ND_m, and dim=sdim=d>1
1728  // - J_{i*}, J's i-th row, is in ( Q_{k-1,k,k}, Q_{k,k-1,k}, Q_{k,k,k-1} )
1729  // - adj(J)_{*j} is in ( Q_{s,s-1,s-1}, Q_{s-1,s,s-1}, Q_{s-1,s-1,s} )
1730  // with s = (d-1)*k
1731  // - adj(J)^T grad_hat(v_hat) is in Q_{(d-1)*k+(l-1)}
1732  // - u_hat is in ( Q_{m-1,m,m}, Q_{m,m-1,m}, Q_{m,m,m-1} )
1733  // - adj(J)^T u_hat is in Q_{(d-1)*k+(m-1)}
1734  // - and u_hat^T adj(J) adj(J)^T grad_hat(v_hat) is in Q_n with
1735  // n = 2*(d-1)*k+(l-1)+(m-1)
1736  //
1737  // In the next formula we use the expressions for n with k=1, which means
1738  // that the term Q/det(J) is disregarded:
1739  int ir_order = (trial_fe.Space() == FunctionSpace::Pk) ?
1740  (trial_fe.GetOrder() + test_fe.GetOrder() - 1) :
1741  (trial_fe.GetOrder() + test_fe.GetOrder() + 2*(dim-2));
1742  ir = &IntRules.Get(trial_fe.GetGeomType(), ir_order);
1743  }
1744 
1745  elmat = 0.0;
1746  for (i = 0; i < ir->GetNPoints(); i++)
1747  {
1748  const IntegrationPoint &ip = ir->IntPoint(i);
1749  test_fe.CalcDShape(ip, dshape);
1750 
1751  Trans.SetIntPoint(&ip);
1752  CalcAdjugate(Trans.Jacobian(), invdfdx);
1753  Mult(dshape, invdfdx, dshapedxt);
1754 
1755  trial_fe.CalcVShape(Trans, vshape);
1756 
1757  double w = ip.weight;
1758 
1759  if (Q)
1760  {
1761  w *= Q->Eval(Trans, ip);
1762  }
1763  dshapedxt *= -w;
1764 
1765  AddMultABt(dshapedxt, vshape, elmat);
1766  }
1767 }
1768 
1769 void VectorFECurlIntegrator::AssembleElementMatrix2(
1770  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1771  ElementTransformation &Trans, DenseMatrix &elmat)
1772 {
1773  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
1774  int dim = trial_fe.GetDim();
1775  int dimc = (dim == 3) ? 3 : 1;
1776 
1777  MFEM_ASSERT(trial_fe.GetMapType() == mfem::FiniteElement::H_CURL ||
1779  "At least one of the finite elements must be in H(Curl)");
1780 
1781  int curl_nd, vec_nd;
1782  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1783  {
1784  curl_nd = trial_nd;
1785  vec_nd = test_nd;
1786  }
1787  else
1788  {
1789  curl_nd = test_nd;
1790  vec_nd = trial_nd;
1791  }
1792 
1793 #ifdef MFEM_THREAD_SAFE
1794  DenseMatrix curlshapeTrial(curl_nd, dimc);
1795  DenseMatrix curlshapeTrial_dFT(curl_nd, dimc);
1796  DenseMatrix vshapeTest(vec_nd, dimc);
1797 #else
1798  curlshapeTrial.SetSize(curl_nd, dimc);
1799  curlshapeTrial_dFT.SetSize(curl_nd, dimc);
1800  vshapeTest.SetSize(vec_nd, dimc);
1801 #endif
1802  Vector shapeTest(vshapeTest.GetData(), vec_nd);
1803 
1804  elmat.SetSize(test_nd, trial_nd);
1805 
1806  const IntegrationRule *ir = IntRule;
1807  if (ir == NULL)
1808  {
1809  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
1810  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1811  }
1812 
1813  elmat = 0.0;
1814  for (i = 0; i < ir->GetNPoints(); i++)
1815  {
1816  const IntegrationPoint &ip = ir->IntPoint(i);
1817 
1818  Trans.SetIntPoint(&ip);
1819  if (dim == 3)
1820  {
1821  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1822  {
1823  trial_fe.CalcCurlShape(ip, curlshapeTrial);
1824  test_fe.CalcVShape(Trans, vshapeTest);
1825  }
1826  else
1827  {
1828  test_fe.CalcCurlShape(ip, curlshapeTrial);
1829  trial_fe.CalcVShape(Trans, vshapeTest);
1830  }
1831  MultABt(curlshapeTrial, Trans.Jacobian(), curlshapeTrial_dFT);
1832  }
1833  else
1834  {
1835  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1836  {
1837  trial_fe.CalcCurlShape(ip, curlshapeTrial_dFT);
1838  test_fe.CalcShape(ip, shapeTest);
1839  }
1840  else
1841  {
1842  test_fe.CalcCurlShape(ip, curlshapeTrial_dFT);
1843  trial_fe.CalcShape(ip, shapeTest);
1844  }
1845  }
1846 
1847  double w = ip.weight;
1848 
1849  if (Q)
1850  {
1851  w *= Q->Eval(Trans, ip);
1852  }
1853  // Note: shapeTest points to the same data as vshapeTest
1854  vshapeTest *= w;
1855  if ( trial_fe.GetMapType() == mfem::FiniteElement::H_CURL )
1856  {
1857  AddMultABt(vshapeTest, curlshapeTrial_dFT, elmat);
1858  }
1859  else
1860  {
1861  AddMultABt(curlshapeTrial_dFT, vshapeTest, elmat);
1862  }
1863  }
1864 }
1865 
1866 void DerivativeIntegrator::AssembleElementMatrix2 (
1867  const FiniteElement &trial_fe,
1868  const FiniteElement &test_fe,
1869  ElementTransformation &Trans,
1870  DenseMatrix &elmat)
1871 {
1872  int dim = trial_fe.GetDim();
1873  int trial_nd = trial_fe.GetDof();
1874  int test_nd = test_fe.GetDof();
1875  int spaceDim = Trans.GetSpaceDim();
1876 
1877  int i, l;
1878  double det;
1879 
1880  elmat.SetSize (test_nd,trial_nd);
1881  dshape.SetSize (trial_nd,dim);
1882  dshapedxt.SetSize(trial_nd, spaceDim);
1883  dshapedxi.SetSize(trial_nd);
1884  invdfdx.SetSize(dim, spaceDim);
1885  shape.SetSize (test_nd);
1886 
1887  const IntegrationRule *ir = IntRule;
1888  if (ir == NULL)
1889  {
1890  int order;
1891  if (trial_fe.Space() == FunctionSpace::Pk)
1892  {
1893  order = trial_fe.GetOrder() + test_fe.GetOrder() - 1;
1894  }
1895  else
1896  {
1897  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim;
1898  }
1899 
1900  if (trial_fe.Space() == FunctionSpace::rQk)
1901  {
1902  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1903  }
1904  else
1905  {
1906  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1907  }
1908  }
1909 
1910  elmat = 0.0;
1911  for (i = 0; i < ir->GetNPoints(); i++)
1912  {
1913  const IntegrationPoint &ip = ir->IntPoint(i);
1914 
1915  trial_fe.CalcDShape(ip, dshape);
1916 
1917  Trans.SetIntPoint (&ip);
1918  CalcInverse (Trans.Jacobian(), invdfdx);
1919  det = Trans.Weight();
1920  Mult (dshape, invdfdx, dshapedxt);
1921 
1922  test_fe.CalcShape(ip, shape);
1923 
1924  for (l = 0; l < trial_nd; l++)
1925  {
1926  dshapedxi(l) = dshapedxt(l,xi);
1927  }
1928 
1929  shape *= Q->Eval(Trans,ip) * det * ip.weight;
1930  AddMultVWt (shape, dshapedxi, elmat);
1931  }
1932 }
1933 
1934 void CurlCurlIntegrator::AssembleElementMatrix
1936  DenseMatrix &elmat )
1937 {
1938  int nd = el.GetDof();
1939  dim = el.GetDim();
1940  int dimc = el.GetCurlDim();
1941  double w;
1942 
1943 #ifdef MFEM_THREAD_SAFE
1944  Vector D;
1945  DenseMatrix curlshape(nd,dimc), curlshape_dFt(nd,dimc), M;
1946 #else
1947  curlshape.SetSize(nd,dimc);
1948  curlshape_dFt.SetSize(nd,dimc);
1949 #endif
1950  elmat.SetSize(nd);
1951  if (MQ) { M.SetSize(dimc); }
1952  if (DQ) { D.SetSize(dimc); }
1953 
1954  const IntegrationRule *ir = IntRule;
1955  if (ir == NULL)
1956  {
1957  int order;
1958  if (el.Space() == FunctionSpace::Pk)
1959  {
1960  order = 2*el.GetOrder() - 2;
1961  }
1962  else
1963  {
1964  order = 2*el.GetOrder();
1965  }
1966 
1967  ir = &IntRules.Get(el.GetGeomType(), order);
1968  }
1969 
1970  elmat = 0.0;
1971  for (int i = 0; i < ir->GetNPoints(); i++)
1972  {
1973  const IntegrationPoint &ip = ir->IntPoint(i);
1974 
1975  Trans.SetIntPoint (&ip);
1976 
1977  w = ip.weight * Trans.Weight();
1978  el.CalcPhysCurlShape(Trans, curlshape_dFt);
1979 
1980  if (MQ)
1981  {
1982  MQ->Eval(M, Trans, ip);
1983  M *= w;
1984  Mult(curlshape_dFt, M, curlshape);
1985  AddMultABt(curlshape, curlshape_dFt, elmat);
1986  }
1987  else if (DQ)
1988  {
1989  DQ->Eval(D, Trans, ip);
1990  D *= w;
1991  AddMultADAt(curlshape_dFt, D, elmat);
1992  }
1993  else if (Q)
1994  {
1995  w *= Q->Eval(Trans, ip);
1996  AddMult_a_AAt(w, curlshape_dFt, elmat);
1997  }
1998  else
1999  {
2000  AddMult_a_AAt(w, curlshape_dFt, elmat);
2001  }
2002  }
2003 }
2004 
2005 void CurlCurlIntegrator
2006 ::ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans,
2007  Vector &u, const FiniteElement &fluxelem, Vector &flux,
2008  bool with_coef)
2009 {
2010 #ifdef MFEM_THREAD_SAFE
2011  DenseMatrix projcurl;
2012 #endif
2013 
2014  fluxelem.ProjectCurl(el, Trans, projcurl);
2015 
2016  flux.SetSize(projcurl.Height());
2017  projcurl.Mult(u, flux);
2018 
2019  // TODO: Q, wcoef?
2020 }
2021 
2022 double CurlCurlIntegrator::ComputeFluxEnergy(const FiniteElement &fluxelem,
2023  ElementTransformation &Trans,
2024  Vector &flux, Vector *d_energy)
2025 {
2026  int nd = fluxelem.GetDof();
2027  dim = fluxelem.GetDim();
2028 
2029 #ifdef MFEM_THREAD_SAFE
2030  DenseMatrix vshape;
2031 #endif
2032  vshape.SetSize(nd, dim);
2033  pointflux.SetSize(dim);
2034  if (d_energy) { vec.SetSize(dim); }
2035 
2036  int order = 2 * fluxelem.GetOrder(); // <--
2037  const IntegrationRule &ir = IntRules.Get(fluxelem.GetGeomType(), order);
2038 
2039  double energy = 0.0;
2040  if (d_energy) { *d_energy = 0.0; }
2041 
2042  Vector* pfluxes = NULL;
2043  if (d_energy)
2044  {
2045  pfluxes = new Vector[ir.GetNPoints()];
2046  }
2047 
2048  for (int i = 0; i < ir.GetNPoints(); i++)
2049  {
2050  const IntegrationPoint &ip = ir.IntPoint(i);
2051  Trans.SetIntPoint(&ip);
2052 
2053  fluxelem.CalcVShape(Trans, vshape);
2054  // fluxelem.CalcVShape(ip, vshape);
2055  vshape.MultTranspose(flux, pointflux);
2056 
2057  double w = Trans.Weight() * ip.weight;
2058 
2059  double e = w * (pointflux * pointflux);
2060 
2061  if (Q)
2062  {
2063  // TODO
2064  }
2065 
2066  energy += e;
2067 
2068 #if ANISO_EXPERIMENTAL
2069  if (d_energy)
2070  {
2071  pfluxes[i].SetSize(dim);
2072  Trans.Jacobian().MultTranspose(pointflux, pfluxes[i]);
2073 
2074  /*
2075  DenseMatrix Jadj(dim, dim);
2076  CalcAdjugate(Trans.Jacobian(), Jadj);
2077  pfluxes[i].SetSize(dim);
2078  Jadj.Mult(pointflux, pfluxes[i]);
2079  */
2080 
2081  // pfluxes[i] = pointflux;
2082  }
2083 #endif
2084  }
2085 
2086  if (d_energy)
2087  {
2088 #if ANISO_EXPERIMENTAL
2089  *d_energy = 0.0;
2090  Vector tmp;
2091 
2092  int n = (int) round(pow(ir.GetNPoints(), 1.0/3.0));
2093  MFEM_ASSERT(n*n*n == ir.GetNPoints(), "");
2094 
2095  // hack: get total variation of 'pointflux' in the x,y,z directions
2096  for (int k = 0; k < n; k++)
2097  for (int l = 0; l < n; l++)
2098  for (int m = 0; m < n; m++)
2099  {
2100  Vector &vec = pfluxes[(k*n + l)*n + m];
2101  if (m > 0)
2102  {
2103  tmp = vec; tmp -= pfluxes[(k*n + l)*n + (m-1)];
2104  (*d_energy)[0] += (tmp * tmp);
2105  }
2106  if (l > 0)
2107  {
2108  tmp = vec; tmp -= pfluxes[(k*n + (l-1))*n + m];
2109  (*d_energy)[1] += (tmp * tmp);
2110  }
2111  if (k > 0)
2112  {
2113  tmp = vec; tmp -= pfluxes[((k-1)*n + l)*n + m];
2114  (*d_energy)[2] += (tmp * tmp);
2115  }
2116  }
2117 #else
2118  *d_energy = 1.0;
2119 #endif
2120 
2121  delete [] pfluxes;
2122  }
2123 
2124  return energy;
2125 }
2126 
2127 void VectorCurlCurlIntegrator::AssembleElementMatrix(
2128  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
2129 {
2130  int dim = el.GetDim();
2131  int dof = el.GetDof();
2132  int cld = (dim*(dim-1))/2;
2133 
2134 #ifdef MFEM_THREAD_SAFE
2135  DenseMatrix dshape_hat(dof, dim), dshape(dof, dim);
2136  DenseMatrix curlshape(dim*dof, cld), Jadj(dim);
2137 #else
2138  dshape_hat.SetSize(dof, dim);
2139  dshape.SetSize(dof, dim);
2140  curlshape.SetSize(dim*dof, cld);
2141  Jadj.SetSize(dim);
2142 #endif
2143 
2144  const IntegrationRule *ir = IntRule;
2145  if (ir == NULL)
2146  {
2147  // use the same integration rule as diffusion
2148  int order = 2 * Trans.OrderGrad(&el);
2149  ir = &IntRules.Get(el.GetGeomType(), order);
2150  }
2151 
2152  elmat.SetSize(dof*dim);
2153  elmat = 0.0;
2154  for (int i = 0; i < ir->GetNPoints(); i++)
2155  {
2156  const IntegrationPoint &ip = ir->IntPoint(i);
2157  el.CalcDShape(ip, dshape_hat);
2158 
2159  Trans.SetIntPoint(&ip);
2160  CalcAdjugate(Trans.Jacobian(), Jadj);
2161  double w = ip.weight / Trans.Weight();
2162 
2163  Mult(dshape_hat, Jadj, dshape);
2164  dshape.GradToCurl(curlshape);
2165 
2166  if (Q)
2167  {
2168  w *= Q->Eval(Trans, ip);
2169  }
2170 
2171  AddMult_a_AAt(w, curlshape, elmat);
2172  }
2173 }
2174 
2175 double VectorCurlCurlIntegrator::GetElementEnergy(
2176  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
2177 {
2178  int dim = el.GetDim();
2179  int dof = el.GetDof();
2180 
2181 #ifdef MFEM_THREAD_SAFE
2182  DenseMatrix dshape_hat(dof, dim), Jadj(dim), grad_hat(dim), grad(dim);
2183 #else
2184  dshape_hat.SetSize(dof, dim);
2185 
2186  Jadj.SetSize(dim);
2187  grad_hat.SetSize(dim);
2188  grad.SetSize(dim);
2189 #endif
2190  DenseMatrix elfun_mat(elfun.GetData(), dof, dim);
2191 
2192  const IntegrationRule *ir = IntRule;
2193  if (ir == NULL)
2194  {
2195  // use the same integration rule as diffusion
2196  int order = 2 * Tr.OrderGrad(&el);
2197  ir = &IntRules.Get(el.GetGeomType(), order);
2198  }
2199 
2200  double energy = 0.;
2201  for (int i = 0; i < ir->GetNPoints(); i++)
2202  {
2203  const IntegrationPoint &ip = ir->IntPoint(i);
2204  el.CalcDShape(ip, dshape_hat);
2205 
2206  MultAtB(elfun_mat, dshape_hat, grad_hat);
2207 
2208  Tr.SetIntPoint(&ip);
2209  CalcAdjugate(Tr.Jacobian(), Jadj);
2210  double w = ip.weight / Tr.Weight();
2211 
2212  Mult(grad_hat, Jadj, grad);
2213 
2214  if (dim == 2)
2215  {
2216  double curl = grad(0,1) - grad(1,0);
2217  w *= curl * curl;
2218  }
2219  else
2220  {
2221  double curl_x = grad(2,1) - grad(1,2);
2222  double curl_y = grad(0,2) - grad(2,0);
2223  double curl_z = grad(1,0) - grad(0,1);
2224  w *= curl_x * curl_x + curl_y * curl_y + curl_z * curl_z;
2225  }
2226 
2227  if (Q)
2228  {
2229  w *= Q->Eval(Tr, ip);
2230  }
2231 
2232  energy += w;
2233  }
2234 
2235  elfun_mat.ClearExternalData();
2236 
2237  return 0.5 * energy;
2238 }
2239 
2240 
2241 void VectorFEMassIntegrator::AssembleElementMatrix(
2242  const FiniteElement &el,
2243  ElementTransformation &Trans,
2244  DenseMatrix &elmat)
2245 {
2246  int dof = el.GetDof();
2247  int spaceDim = Trans.GetSpaceDim();
2248  int vdim = std::max(spaceDim, el.GetVDim());
2249 
2250  double w;
2251 
2252 #ifdef MFEM_THREAD_SAFE
2253  Vector D(DQ ? DQ->GetVDim() : 0);
2254  DenseMatrix trial_vshape(dof, vdim);
2255  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2256 #else
2257  trial_vshape.SetSize(dof, vdim);
2258  D.SetSize(DQ ? DQ->GetVDim() : 0);
2259  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2260 #endif
2261  DenseMatrix tmp(trial_vshape.Height(), K.Width());
2262 
2263  elmat.SetSize(dof);
2264  elmat = 0.0;
2265 
2266  const IntegrationRule *ir = IntRule;
2267  if (ir == NULL)
2268  {
2269  // int order = 2 * el.GetOrder();
2270  int order = Trans.OrderW() + 2 * el.GetOrder();
2271  ir = &IntRules.Get(el.GetGeomType(), order);
2272  }
2273 
2274  for (int i = 0; i < ir->GetNPoints(); i++)
2275  {
2276  const IntegrationPoint &ip = ir->IntPoint(i);
2277 
2278  Trans.SetIntPoint (&ip);
2279 
2280  el.CalcVShape(Trans, trial_vshape);
2281 
2282  w = ip.weight * Trans.Weight();
2283  if (MQ)
2284  {
2285  MQ->Eval(K, Trans, ip);
2286  K *= w;
2287  Mult(trial_vshape,K,tmp);
2288  AddMultABt(tmp,trial_vshape,elmat);
2289  }
2290  else if (DQ)
2291  {
2292  DQ->Eval(D, Trans, ip);
2293  D *= w;
2294  AddMultADAt(trial_vshape, D, elmat);
2295  }
2296  else
2297  {
2298  if (Q)
2299  {
2300  w *= Q -> Eval (Trans, ip);
2301  }
2302  AddMult_a_AAt (w, trial_vshape, elmat);
2303  }
2304  }
2305 }
2306 
2307 void VectorFEMassIntegrator::AssembleElementMatrix2(
2308  const FiniteElement &trial_fe, const FiniteElement &test_fe,
2309  ElementTransformation &Trans, DenseMatrix &elmat)
2310 {
2311  if (test_fe.GetRangeType() == FiniteElement::SCALAR
2312  && trial_fe.GetRangeType() == FiniteElement::VECTOR)
2313  {
2314  // assume test_fe is scalar FE and trial_fe is vector FE
2315  int spaceDim = Trans.GetSpaceDim();
2316  int vdim = std::max(spaceDim, trial_fe.GetVDim());
2317  int trial_dof = trial_fe.GetDof();
2318  int test_dof = test_fe.GetDof();
2319  double w;
2320 
2321 #ifdef MFEM_THREAD_SAFE
2322  DenseMatrix trial_vshape(trial_dof, spaceDim);
2323  Vector shape(test_dof);
2324  Vector D(DQ ? DQ->GetVDim() : 0);
2325  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2326 #else
2327  trial_vshape.SetSize(trial_dof, spaceDim);
2328  shape.SetSize(test_dof);
2329  D.SetSize(DQ ? DQ->GetVDim() : 0);
2330  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2331 #endif
2332 
2333  elmat.SetSize(vdim*test_dof, trial_dof);
2334 
2335  const IntegrationRule *ir = IntRule;
2336  if (ir == NULL)
2337  {
2338  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
2339  ir = &IntRules.Get(test_fe.GetGeomType(), order);
2340  }
2341 
2342  elmat = 0.0;
2343  for (int i = 0; i < ir->GetNPoints(); i++)
2344  {
2345  const IntegrationPoint &ip = ir->IntPoint(i);
2346 
2347  Trans.SetIntPoint (&ip);
2348 
2349  trial_fe.CalcVShape(Trans, trial_vshape);
2350  test_fe.CalcShape(ip, shape);
2351 
2352  w = ip.weight * Trans.Weight();
2353  if (DQ)
2354  {
2355  DQ->Eval(D, Trans, ip);
2356  D *= w;
2357  for (int d = 0; d < vdim; d++)
2358  {
2359  for (int j = 0; j < test_dof; j++)
2360  {
2361  for (int k = 0; k < trial_dof; k++)
2362  {
2363  elmat(d * test_dof + j, k) +=
2364  shape(j) * D(d) * trial_vshape(k, d);
2365  }
2366  }
2367  }
2368  }
2369  else if (MQ)
2370  {
2371  MQ->Eval(K, Trans, ip);
2372  K *= w;
2373  for (int d = 0; d < vdim; d++)
2374  {
2375  for (int j = 0; j < test_dof; j++)
2376  {
2377  for (int k = 0; k < trial_dof; k++)
2378  {
2379  double Kv = 0.0;
2380  for (int vd = 0; vd < spaceDim; vd++)
2381  {
2382  Kv += K(d, vd) * trial_vshape(k, vd);
2383  }
2384  elmat(d * test_dof + j, k) += shape(j) * Kv;
2385  }
2386  }
2387  }
2388  }
2389  else
2390  {
2391  if (Q)
2392  {
2393  w *= Q->Eval(Trans, ip);
2394  }
2395  for (int d = 0; d < vdim; d++)
2396  {
2397  for (int j = 0; j < test_dof; j++)
2398  {
2399  for (int k = 0; k < trial_dof; k++)
2400  {
2401  elmat(d * test_dof + j, k) +=
2402  w * shape(j) * trial_vshape(k, d);
2403  }
2404  }
2405  }
2406  }
2407  }
2408  }
2409  else if (test_fe.GetRangeType() == FiniteElement::VECTOR
2410  && trial_fe.GetRangeType() == FiniteElement::VECTOR)
2411  {
2412  // assume both test_fe and trial_fe are vector FE
2413  int spaceDim = Trans.GetSpaceDim();
2414  int trial_vdim = std::max(spaceDim, trial_fe.GetVDim());
2415  int test_vdim = std::max(spaceDim, test_fe.GetVDim());
2416  int trial_dof = trial_fe.GetDof();
2417  int test_dof = test_fe.GetDof();
2418  double w;
2419 
2420 #ifdef MFEM_THREAD_SAFE
2421  DenseMatrix trial_vshape(trial_dof,trial_vdim);
2422  DenseMatrix test_vshape(test_dof,test_vdim);
2423  Vector D(DQ ? DQ->GetVDim() : 0);
2424  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2425 #else
2426  trial_vshape.SetSize(trial_dof,trial_vdim);
2427  test_vshape.SetSize(test_dof,test_vdim);
2428  D.SetSize(DQ ? DQ->GetVDim() : 0);
2429  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
2430 #endif
2431  DenseMatrix tmp(test_vshape.Height(), K.Width());
2432 
2433  elmat.SetSize (test_dof, trial_dof);
2434 
2435  const IntegrationRule *ir = IntRule;
2436  if (ir == NULL)
2437  {
2438  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
2439  ir = &IntRules.Get(test_fe.GetGeomType(), order);
2440  }
2441 
2442  elmat = 0.0;
2443  for (int i = 0; i < ir->GetNPoints(); i++)
2444  {
2445  const IntegrationPoint &ip = ir->IntPoint(i);
2446 
2447  Trans.SetIntPoint (&ip);
2448 
2449  trial_fe.CalcVShape(Trans, trial_vshape);
2450  test_fe.CalcVShape(Trans, test_vshape);
2451 
2452  w = ip.weight * Trans.Weight();
2453  if (MQ)
2454  {
2455  MQ->Eval(K, Trans, ip);
2456  K *= w;
2457  Mult(test_vshape,K,tmp);
2458  AddMultABt(tmp,trial_vshape,elmat);
2459  }
2460  else if (DQ)
2461  {
2462  DQ->Eval(D, Trans, ip);
2463  D *= w;
2464  AddMultADBt(test_vshape,D,trial_vshape,elmat);
2465  }
2466  else
2467  {
2468  if (Q)
2469  {
2470  w *= Q -> Eval (Trans, ip);
2471  }
2472  AddMult_a_ABt(w,test_vshape,trial_vshape,elmat);
2473  }
2474  }
2475  }
2476  else
2477  {
2478  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
2479  " is not implemented for given trial and test bases.");
2480  }
2481 }
2482 
2483 void VectorDivergenceIntegrator::AssembleElementMatrix2(
2484  const FiniteElement &trial_fe,
2485  const FiniteElement &test_fe,
2486  ElementTransformation &Trans,
2487  DenseMatrix &elmat)
2488 {
2489  dim = trial_fe.GetDim();
2490  int trial_dof = trial_fe.GetDof();
2491  int test_dof = test_fe.GetDof();
2492  double c;
2493 
2494  dshape.SetSize (trial_dof, dim);
2495  gshape.SetSize (trial_dof, dim);
2496  Jadj.SetSize (dim);
2497  divshape.SetSize (dim*trial_dof);
2498  shape.SetSize (test_dof);
2499 
2500  elmat.SetSize (test_dof, dim*trial_dof);
2501 
2502  const IntegrationRule *ir = IntRule ? IntRule : &GetRule(trial_fe, test_fe,
2503  Trans);
2504 
2505  elmat = 0.0;
2506 
2507  for (int i = 0; i < ir -> GetNPoints(); i++)
2508  {
2509  const IntegrationPoint &ip = ir->IntPoint(i);
2510 
2511  trial_fe.CalcDShape (ip, dshape);
2512  test_fe.CalcShape (ip, shape);
2513 
2514  Trans.SetIntPoint (&ip);
2515  CalcAdjugate(Trans.Jacobian(), Jadj);
2516 
2517  Mult (dshape, Jadj, gshape);
2518 
2519  gshape.GradToDiv (divshape);
2520 
2521  c = ip.weight;
2522  if (Q)
2523  {
2524  c *= Q -> Eval (Trans, ip);
2525  }
2526 
2527  // elmat += c * shape * divshape ^ t
2528  shape *= c;
2529  AddMultVWt (shape, divshape, elmat);
2530  }
2531 }
2532 
2533 const IntegrationRule &VectorDivergenceIntegrator::GetRule(
2534  const FiniteElement &trial_fe,
2535  const FiniteElement &test_fe,
2536  ElementTransformation &Trans)
2537 {
2538  int order = Trans.OrderGrad(&trial_fe) + test_fe.GetOrder() + Trans.OrderJ();
2539  return IntRules.Get(trial_fe.GetGeomType(), order);
2540 }
2541 
2542 
2543 void DivDivIntegrator::AssembleElementMatrix(
2544  const FiniteElement &el,
2545  ElementTransformation &Trans,
2546  DenseMatrix &elmat)
2547 {
2548  int dof = el.GetDof();
2549  double c;
2550 
2551 #ifdef MFEM_THREAD_SAFE
2552  Vector divshape(dof);
2553 #else
2554  divshape.SetSize(dof);
2555 #endif
2556  elmat.SetSize(dof);
2557 
2558  const IntegrationRule *ir = IntRule;
2559  if (ir == NULL)
2560  {
2561  int order = 2 * el.GetOrder() - 2; // <--- OK for RTk
2562  ir = &IntRules.Get(el.GetGeomType(), order);
2563  }
2564 
2565  elmat = 0.0;
2566 
2567  for (int i = 0; i < ir -> GetNPoints(); i++)
2568  {
2569  const IntegrationPoint &ip = ir->IntPoint(i);
2570 
2571  el.CalcDivShape (ip, divshape);
2572 
2573  Trans.SetIntPoint (&ip);
2574  c = ip.weight / Trans.Weight();
2575 
2576  if (Q)
2577  {
2578  c *= Q -> Eval (Trans, ip);
2579  }
2580 
2581  // elmat += c * divshape * divshape ^ t
2582  AddMult_a_VVt (c, divshape, elmat);
2583  }
2584 }
2585 
2586 
2587 void VectorDiffusionIntegrator::AssembleElementMatrix(
2588  const FiniteElement &el,
2589  ElementTransformation &Trans,
2590  DenseMatrix &elmat)
2591 {
2592  const int dof = el.GetDof();
2593  dim = el.GetDim();
2594  sdim = Trans.GetSpaceDim();
2595 
2596  // If vdim is not set, set it to the space dimension;
2597  vdim = (vdim <= 0) ? sdim : vdim;
2598  const bool square = (dim == sdim);
2599 
2600  if (VQ)
2601  {
2602  vcoeff.SetSize(vdim);
2603  }
2604  else if (MQ)
2605  {
2606  mcoeff.SetSize(vdim);
2607  }
2608 
2609  dshape.SetSize(dof, dim);
2610  dshapedxt.SetSize(dof, sdim);
2611 
2612  elmat.SetSize(vdim * dof);
2613  pelmat.SetSize(dof);
2614 
2615  const IntegrationRule *ir = IntRule;
2616  if (ir == NULL)
2617  {
2618  ir = &DiffusionIntegrator::GetRule(el,el);
2619  }
2620 
2621  elmat = 0.0;
2622 
2623  for (int i = 0; i < ir -> GetNPoints(); i++)
2624  {
2625 
2626  const IntegrationPoint &ip = ir->IntPoint(i);
2627  el.CalcDShape(ip, dshape);
2628 
2629  Trans.SetIntPoint(&ip);
2630  double w = Trans.Weight();
2631  w = ip.weight / (square ? w : w*w*w);
2632  // AdjugateJacobian = / adj(J), if J is square
2633  // \ adj(J^t.J).J^t, otherwise
2634  Mult(dshape, Trans.AdjugateJacobian(), dshapedxt);
2635 
2636  if (VQ)
2637  {
2638  VQ->Eval(vcoeff, Trans, ip);
2639  for (int k = 0; k < vdim; ++k)
2640  {
2641  Mult_a_AAt(w*vcoeff(k), dshapedxt, pelmat);
2642  elmat.AddMatrix(pelmat, dof*k, dof*k);
2643  }
2644  }
2645  else if (MQ)
2646  {
2647  MQ->Eval(mcoeff, Trans, ip);
2648  for (int ii = 0; ii < vdim; ++ii)
2649  {
2650  for (int jj = 0; jj < vdim; ++jj)
2651  {
2652  Mult_a_AAt(w*mcoeff(ii,jj), dshapedxt, pelmat);
2653  elmat.AddMatrix(pelmat, dof*ii, dof*jj);
2654  }
2655  }
2656  }
2657  else
2658  {
2659  if (Q) { w *= Q->Eval(Trans, ip); }
2660  Mult_a_AAt(w, dshapedxt, pelmat);
2661  for (int k = 0; k < vdim; ++k)
2662  {
2663  elmat.AddMatrix(pelmat, dof*k, dof*k);
2664  }
2665  }
2666  }
2667 }
2668 
2669 void VectorDiffusionIntegrator::AssembleElementVector(
2670  const FiniteElement &el, ElementTransformation &Tr,
2671  const Vector &elfun, Vector &elvect)
2672 {
2673  const int dof = el.GetDof();
2674  dim = el.GetDim();
2675  sdim = Tr.GetSpaceDim();
2676 
2677  // If vdim is not set, set it to the space dimension;
2678  vdim = (vdim <= 0) ? sdim : vdim;
2679  const bool square = (dim == sdim);
2680 
2681  if (VQ)
2682  {
2683  vcoeff.SetSize(vdim);
2684  }
2685  else if (MQ)
2686  {
2687  mcoeff.SetSize(vdim);
2688  }
2689 
2690  dshape.SetSize(dof, dim);
2691  dshapedxt.SetSize(dof, dim);
2692  // pelmat.SetSize(dim);
2693 
2694  elvect.SetSize(dim*dof);
2695 
2696  // NOTE: DenseMatrix is in column-major order. This is consistent with
2697  // vectors ordered byNODES. In the resulting DenseMatrix, each column
2698  // corresponds to a particular vdim.
2699  DenseMatrix mat_in(elfun.GetData(), dof, dim);
2700  DenseMatrix mat_out(elvect.GetData(), dof, dim);
2701 
2702  const IntegrationRule *ir = IntRule;
2703  if (ir == NULL)
2704  {
2705  ir = &DiffusionIntegrator::GetRule(el,el);
2706  }
2707 
2708  elvect = 0.0;
2709  for (int i = 0; i < ir->GetNPoints(); i++)
2710  {
2711  const IntegrationPoint &ip = ir->IntPoint(i);
2712  el.CalcDShape(ip, dshape);
2713 
2714  Tr.SetIntPoint(&ip);
2715  double w = Tr.Weight();
2716  w = ip.weight / (square ? w : w*w*w);
2717  Mult(dshape, Tr.AdjugateJacobian(), dshapedxt);
2718  MultAAt(dshapedxt, pelmat);
2719 
2720  if (VQ)
2721  {
2722  VQ->Eval(vcoeff, Tr, ip);
2723  for (int k = 0; k < vdim; ++k)
2724  {
2725  pelmat *= w*vcoeff(k);
2726  const Vector vec_in(mat_in.GetColumn(k), dof);
2727  Vector vec_out(mat_out.GetColumn(k), dof);
2728  pelmat.AddMult(vec_in, vec_out);
2729  }
2730  }
2731  else if (MQ)
2732  {
2733  MQ->Eval(mcoeff, Tr, ip);
2734  for (int ii = 0; ii < vdim; ++ii)
2735  {
2736  Vector vec_out(mat_out.GetColumn(ii), dof);
2737  for (int jj = 0; jj < vdim; ++jj)
2738  {
2739  pelmat *= w*mcoeff(ii,jj);
2740  const Vector vec_in(mat_in.GetColumn(jj), dof);
2741  pelmat.Mult(vec_in, vec_out);
2742  }
2743  }
2744  }
2745  else
2746  {
2747  if (Q) { w *= Q->Eval(Tr, ip); }
2748  pelmat *= w;
2749  for (int k = 0; k < vdim; ++k)
2750  {
2751  const Vector vec_in(mat_in.GetColumn(k), dof);
2752  Vector vec_out(mat_out.GetColumn(k), dof);
2753  pelmat.AddMult(vec_in, vec_out);
2754  }
2755  }
2756  }
2757 }
2758 
2759 
2760 void ElasticityIntegrator::AssembleElementMatrix(
2761  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
2762 {
2763  int dof = el.GetDof();
2764  int dim = el.GetDim();
2765  double w, L, M;
2766 
2767  MFEM_ASSERT(dim == Trans.GetSpaceDim(), "");
2768 
2769 #ifdef MFEM_THREAD_SAFE
2770  DenseMatrix dshape(dof, dim), gshape(dof, dim), pelmat(dof);
2771  Vector divshape(dim*dof);
2772 #else
2773  dshape.SetSize(dof, dim);
2774  gshape.SetSize(dof, dim);
2775  pelmat.SetSize(dof);
2776  divshape.SetSize(dim*dof);
2777 #endif
2778 
2779  elmat.SetSize(dof * dim);
2780 
2781  const IntegrationRule *ir = IntRule;
2782  if (ir == NULL)
2783  {
2784  int order = 2 * Trans.OrderGrad(&el); // correct order?
2785  ir = &IntRules.Get(el.GetGeomType(), order);
2786  }
2787 
2788  elmat = 0.0;
2789 
2790  for (int i = 0; i < ir -> GetNPoints(); i++)
2791  {
2792  const IntegrationPoint &ip = ir->IntPoint(i);
2793 
2794  el.CalcDShape(ip, dshape);
2795 
2796  Trans.SetIntPoint(&ip);
2797  w = ip.weight * Trans.Weight();
2798  Mult(dshape, Trans.InverseJacobian(), gshape);
2799  MultAAt(gshape, pelmat);
2800  gshape.GradToDiv (divshape);
2801 
2802  M = mu->Eval(Trans, ip);
2803  if (lambda)
2804  {
2805  L = lambda->Eval(Trans, ip);
2806  }
2807  else
2808  {
2809  L = q_lambda * M;
2810  M = q_mu * M;
2811  }
2812 
2813  if (L != 0.0)
2814  {
2815  AddMult_a_VVt(L * w, divshape, elmat);
2816  }
2817 
2818  if (M != 0.0)
2819  {
2820  for (int d = 0; d < dim; d++)
2821  {
2822  for (int k = 0; k < dof; k++)
2823  for (int l = 0; l < dof; l++)
2824  {
2825  elmat (dof*d+k, dof*d+l) += (M * w) * pelmat(k, l);
2826  }
2827  }
2828  for (int ii = 0; ii < dim; ii++)
2829  for (int jj = 0; jj < dim; jj++)
2830  {
2831  for (int kk = 0; kk < dof; kk++)
2832  for (int ll = 0; ll < dof; ll++)
2833  {
2834  elmat(dof*ii+kk, dof*jj+ll) +=
2835  (M * w) * gshape(kk, jj) * gshape(ll, ii);
2836  }
2837  }
2838  }
2839  }
2840 }
2841 
2842 void ElasticityIntegrator::ComputeElementFlux(
2843  const mfem::FiniteElement &el, ElementTransformation &Trans,
2844  Vector &u, const mfem::FiniteElement &fluxelem, Vector &flux,
2845  bool with_coef)
2846 {
2847  const int dof = el.GetDof();
2848  const int dim = el.GetDim();
2849  const int tdim = dim*(dim+1)/2; // num. entries in a symmetric tensor
2850  double L, M;
2851 
2852  MFEM_ASSERT(dim == 2 || dim == 3,
2853  "dimension is not supported: dim = " << dim);
2854  MFEM_ASSERT(dim == Trans.GetSpaceDim(), "");
2855  MFEM_ASSERT(fluxelem.GetMapType() == FiniteElement::VALUE, "");
2856  MFEM_ASSERT(dynamic_cast<const NodalFiniteElement*>(&fluxelem), "");
2857 
2858 #ifdef MFEM_THREAD_SAFE
2859  DenseMatrix dshape(dof, dim);
2860 #else
2861  dshape.SetSize(dof, dim);
2862 #endif
2863 
2864  double gh_data[9], grad_data[9];
2865  DenseMatrix gh(gh_data, dim, dim);
2866  DenseMatrix grad(grad_data, dim, dim);
2867 
2868  const IntegrationRule &ir = fluxelem.GetNodes();
2869  const int fnd = ir.GetNPoints();
2870  flux.SetSize(fnd * tdim);
2871 
2872  DenseMatrix loc_data_mat(u.GetData(), dof, dim);
2873  for (int i = 0; i < fnd; i++)
2874  {
2875  const IntegrationPoint &ip = ir.IntPoint(i);
2876  el.CalcDShape(ip, dshape);
2877  MultAtB(loc_data_mat, dshape, gh);
2878 
2879  Trans.SetIntPoint(&ip);
2880  Mult(gh, Trans.InverseJacobian(), grad);
2881 
2882  M = mu->Eval(Trans, ip);
2883  if (lambda)
2884  {
2885  L = lambda->Eval(Trans, ip);
2886  }
2887  else
2888  {
2889  L = q_lambda * M;
2890  M = q_mu * M;
2891  }
2892 
2893  // stress = 2*M*e(u) + L*tr(e(u))*I, where
2894  // e(u) = (1/2)*(grad(u) + grad(u)^T)
2895  const double M2 = 2.0*M;
2896  if (dim == 2)
2897  {
2898  L *= (grad(0,0) + grad(1,1));
2899  // order of the stress entries: s_xx, s_yy, s_xy
2900  flux(i+fnd*0) = M2*grad(0,0) + L;
2901  flux(i+fnd*1) = M2*grad(1,1) + L;
2902  flux(i+fnd*2) = M*(grad(0,1) + grad(1,0));
2903  }
2904  else if (dim == 3)
2905  {
2906  L *= (grad(0,0) + grad(1,1) + grad(2,2));
2907  // order of the stress entries: s_xx, s_yy, s_zz, s_xy, s_xz, s_yz
2908  flux(i+fnd*0) = M2*grad(0,0) + L;
2909  flux(i+fnd*1) = M2*grad(1,1) + L;
2910  flux(i+fnd*2) = M2*grad(2,2) + L;
2911  flux(i+fnd*3) = M*(grad(0,1) + grad(1,0));
2912  flux(i+fnd*4) = M*(grad(0,2) + grad(2,0));
2913  flux(i+fnd*5) = M*(grad(1,2) + grad(2,1));
2914  }
2915  }
2916 }
2917 
2918 double ElasticityIntegrator::ComputeFluxEnergy(const FiniteElement &fluxelem,
2919  ElementTransformation &Trans,
2920  Vector &flux, Vector *d_energy)
2921 {
2922  const int dof = fluxelem.GetDof();
2923  const int dim = fluxelem.GetDim();
2924  const int tdim = dim*(dim+1)/2; // num. entries in a symmetric tensor
2925  double L, M;
2926 
2927  // The MFEM_ASSERT constraints in ElasticityIntegrator::ComputeElementFlux
2928  // are assumed here too.
2929  MFEM_ASSERT(d_energy == NULL, "anisotropic estimates are not supported");
2930  MFEM_ASSERT(flux.Size() == dof*tdim, "invalid 'flux' vector");
2931 
2932 #ifndef MFEM_THREAD_SAFE
2933  shape.SetSize(dof);
2934 #else
2935  Vector shape(dof);
2936 #endif
2937  double pointstress_data[6];
2938  Vector pointstress(pointstress_data, tdim);
2939 
2940  // View of the 'flux' vector as a (dof x tdim) matrix
2941  DenseMatrix flux_mat(flux.GetData(), dof, tdim);
2942 
2943  // Use the same integration rule as in AssembleElementMatrix, replacing 'el'
2944  // with 'fluxelem' when 'IntRule' is not set.
2945  // Should we be using a different (more accurate) rule here?
2946  const IntegrationRule *ir = IntRule;
2947  if (ir == NULL)
2948  {
2949  int order = 2 * Trans.OrderGrad(&fluxelem);
2950  ir = &IntRules.Get(fluxelem.GetGeomType(), order);
2951  }
2952 
2953  double energy = 0.0;
2954 
2955  for (int i = 0; i < ir->GetNPoints(); i++)
2956  {
2957  const IntegrationPoint &ip = ir->IntPoint(i);
2958  fluxelem.CalcShape(ip, shape);
2959 
2960  flux_mat.MultTranspose(shape, pointstress);
2961 
2962  Trans.SetIntPoint(&ip);
2963  double w = Trans.Weight() * ip.weight;
2964 
2965  M = mu->Eval(Trans, ip);
2966  if (lambda)
2967  {
2968  L = lambda->Eval(Trans, ip);
2969  }
2970  else
2971  {
2972  L = q_lambda * M;
2973  M = q_mu * M;
2974  }
2975 
2976  // The strain energy density at a point is given by (1/2)*(s : e) where s
2977  // and e are the stress and strain tensors, respectively. Since we only
2978  // have the stress, we need to compute the strain from the stress:
2979  // s = 2*mu*e + lambda*tr(e)*I
2980  // Taking trace on both sides we find:
2981  // tr(s) = 2*mu*tr(e) + lambda*tr(e)*dim = (2*mu + dim*lambda)*tr(e)
2982  // which gives:
2983  // tr(e) = tr(s)/(2*mu + dim*lambda)
2984  // Then from the first identity above we can find the strain:
2985  // e = (1/(2*mu))*(s - lambda*tr(e)*I)
2986 
2987  double pt_e; // point strain energy density
2988  const double *s = pointstress_data;
2989  if (dim == 2)
2990  {
2991  // s entries: s_xx, s_yy, s_xy
2992  const double tr_e = (s[0] + s[1])/(2*(M + L));
2993  L *= tr_e;
2994  pt_e = (0.25/M)*(s[0]*(s[0] - L) + s[1]*(s[1] - L) + 2*s[2]*s[2]);
2995  }
2996  else // (dim == 3)
2997  {
2998  // s entries: s_xx, s_yy, s_zz, s_xy, s_xz, s_yz
2999  const double tr_e = (s[0] + s[1] + s[2])/(2*M + 3*L);
3000  L *= tr_e;
3001  pt_e = (0.25/M)*(s[0]*(s[0] - L) + s[1]*(s[1] - L) + s[2]*(s[2] - L) +
3002  2*(s[3]*s[3] + s[4]*s[4] + s[5]*s[5]));
3003  }
3004 
3005  energy += w * pt_e;
3006  }
3007 
3008  return energy;
3009 }
3010 
3011 void DGTraceIntegrator::AssembleFaceMatrix(const FiniteElement &el1,
3012  const FiniteElement &el2,
3014  DenseMatrix &elmat)
3015 {
3016  int ndof1, ndof2;
3017 
3018  double un, a, b, w;
3019 
3020  dim = el1.GetDim();
3021  ndof1 = el1.GetDof();
3022  Vector vu(dim), nor(dim);
3023 
3024  if (Trans.Elem2No >= 0)
3025  {
3026  ndof2 = el2.GetDof();
3027  }
3028  else
3029  {
3030  ndof2 = 0;
3031  }
3032 
3033  shape1.SetSize(ndof1);
3034  shape2.SetSize(ndof2);
3035  elmat.SetSize(ndof1 + ndof2);
3036  elmat = 0.0;
3037 
3038  const IntegrationRule *ir = IntRule;
3039  if (ir == NULL)
3040  {
3041  int order;
3042  // Assuming order(u)==order(mesh)
3043  if (Trans.Elem2No >= 0)
3044  order = (min(Trans.Elem1->OrderW(), Trans.Elem2->OrderW()) +
3045  2*max(el1.GetOrder(), el2.GetOrder()));
3046  else
3047  {
3048  order = Trans.Elem1->OrderW() + 2*el1.GetOrder();
3049  }
3050  if (el1.Space() == FunctionSpace::Pk)
3051  {
3052  order++;
3053  }
3054  ir = &IntRules.Get(Trans.GetGeometryType(), order);
3055  }
3056 
3057  for (int p = 0; p < ir->GetNPoints(); p++)
3058  {
3059  const IntegrationPoint &ip = ir->IntPoint(p);
3060 
3061  // Set the integration point in the face and the neighboring elements
3062  Trans.SetAllIntPoints(&ip);
3063 
3064  // Access the neighboring elements' integration points
3065  // Note: eip2 will only contain valid data if Elem2 exists
3066  const IntegrationPoint &eip1 = Trans.GetElement1IntPoint();
3067  const IntegrationPoint &eip2 = Trans.GetElement2IntPoint();
3068 
3069  el1.CalcShape(eip1, shape1);
3070 
3071  u->Eval(vu, *Trans.Elem1, eip1);
3072 
3073  if (dim == 1)
3074  {
3075  nor(0) = 2*eip1.x - 1.0;
3076  }
3077  else
3078  {
3079  CalcOrtho(Trans.Jacobian(), nor);
3080  }
3081 
3082  un = vu * nor;
3083  a = 0.5 * alpha * un;
3084  b = beta * fabs(un);
3085  // note: if |alpha/2|==|beta| then |a|==|b|, i.e. (a==b) or (a==-b)
3086  // and therefore two blocks in the element matrix contribution
3087  // (from the current quadrature point) are 0
3088 
3089  if (rho)
3090  {
3091  double rho_p;
3092  if (un >= 0.0 && ndof2)
3093  {
3094  rho_p = rho->Eval(*Trans.Elem2, eip2);
3095  }
3096  else
3097  {
3098  rho_p = rho->Eval(*Trans.Elem1, eip1);
3099  }
3100  a *= rho_p;
3101  b *= rho_p;
3102  }
3103 
3104  w = ip.weight * (a+b);
3105  if (w != 0.0)
3106  {
3107  for (int i = 0; i < ndof1; i++)
3108  for (int j = 0; j < ndof1; j++)
3109  {
3110  elmat(i, j) += w * shape1(i) * shape1(j);
3111  }
3112  }
3113 
3114  if (ndof2)
3115  {
3116  el2.CalcShape(eip2, shape2);
3117 
3118  if (w != 0.0)
3119  for (int i = 0; i < ndof2; i++)
3120  for (int j = 0; j < ndof1; j++)
3121  {
3122  elmat(ndof1+i, j) -= w * shape2(i) * shape1(j);
3123  }
3124 
3125  w = ip.weight * (b-a);
3126  if (w != 0.0)
3127  {
3128  for (int i = 0; i < ndof2; i++)
3129  for (int j = 0; j < ndof2; j++)
3130  {
3131  elmat(ndof1+i, ndof1+j) += w * shape2(i) * shape2(j);
3132  }
3133 
3134  for (int i = 0; i < ndof1; i++)
3135  for (int j = 0; j < ndof2; j++)
3136  {
3137  elmat(i, ndof1+j) -= w * shape1(i) * shape2(j);
3138  }
3139  }
3140  }
3141  }
3142 }
3143 
3144 
3145 const IntegrationRule &DGTraceIntegrator::GetRule(
3146  Geometry::Type geom, int order, FaceElementTransformations &T)
3147 {
3148  int int_order = T.Elem1->OrderW() + 2*order;
3149  return IntRules.Get(geom, int_order);
3150 }
3151 
3152 void DGDiffusionIntegrator::AssembleFaceMatrix(
3153  const FiniteElement &el1, const FiniteElement &el2,
3154  FaceElementTransformations &Trans, DenseMatrix &elmat)
3155 {
3156  int dim, ndof1, ndof2, ndofs;
3157  bool kappa_is_nonzero = (kappa != 0.);
3158  double w, wq = 0.0;
3159 
3160  dim = el1.GetDim();
3161  ndof1 = el1.GetDof();
3162 
3163  nor.SetSize(dim);
3164  nh.SetSize(dim);
3165  ni.SetSize(dim);
3166  adjJ.SetSize(dim);
3167  if (MQ)
3168  {
3169  mq.SetSize(dim);
3170  }
3171 
3172  shape1.SetSize(ndof1);
3173  dshape1.SetSize(ndof1, dim);
3174  dshape1dn.SetSize(ndof1);
3175  if (Trans.Elem2No >= 0)
3176  {
3177  ndof2 = el2.GetDof();
3178  shape2.SetSize(ndof2);
3179  dshape2.SetSize(ndof2, dim);
3180  dshape2dn.SetSize(ndof2);
3181  }
3182  else
3183  {
3184  ndof2 = 0;
3185  }
3186 
3187  ndofs = ndof1 + ndof2;
3188  elmat.SetSize(ndofs);
3189  elmat = 0.0;
3190  if (kappa_is_nonzero)
3191  {
3192  jmat.SetSize(ndofs);
3193  jmat = 0.;
3194  }
3195 
3196  const IntegrationRule *ir = IntRule;
3197  if (ir == NULL)
3198  {
3199  // a simple choice for the integration order; is this OK?
3200  int order;
3201  if (ndof2)
3202  {
3203  order = 2*max(el1.GetOrder(), el2.GetOrder());
3204  }
3205  else
3206  {
3207  order = 2*el1.GetOrder();
3208  }
3209  ir = &IntRules.Get(Trans.GetGeometryType(), order);
3210  }
3211 
3212  // assemble: < {(Q \nabla u).n},[v] > --> elmat
3213  // kappa < {h^{-1} Q} [u],[v] > --> jmat
3214  for (int p = 0; p < ir->GetNPoints(); p++)
3215  {
3216  const IntegrationPoint &ip = ir->IntPoint(p);
3217 
3218  // Set the integration point in the face and the neighboring elements
3219  Trans.SetAllIntPoints(&ip);
3220 
3221  // Access the neighboring elements' integration points
3222  // Note: eip2 will only contain valid data if Elem2 exists
3223  const IntegrationPoint &eip1 = Trans.GetElement1IntPoint();
3224  const IntegrationPoint &eip2 = Trans.GetElement2IntPoint();
3225 
3226  if (dim == 1)
3227  {
3228  nor(0) = 2*eip1.x - 1.0;
3229  }
3230  else
3231  {
3232  CalcOrtho(Trans.Jacobian(), nor);
3233  }
3234 
3235  el1.CalcShape(eip1, shape1);
3236  el1.CalcDShape(eip1, dshape1);
3237  w = ip.weight/Trans.Elem1->Weight();
3238  if (ndof2)
3239  {
3240  w /= 2;
3241  }
3242  if (!MQ)
3243  {
3244  if (Q)
3245  {
3246  w *= Q->Eval(*Trans.Elem1, eip1);
3247  }
3248  ni.Set(w, nor);
3249  }
3250  else
3251  {
3252  nh.Set(w, nor);
3253  MQ->Eval(mq, *Trans.Elem1, eip1);
3254  mq.MultTranspose(nh, ni);
3255  }
3256  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
3257  adjJ.Mult(ni, nh);
3258  if (kappa_is_nonzero)
3259  {
3260  wq = ni * nor;
3261  }
3262  // Note: in the jump term, we use 1/h1 = |nor|/det(J1) which is
3263  // independent of Loc1 and always gives the size of element 1 in
3264  // direction perpendicular to the face. Indeed, for linear transformation
3265  // |nor|=measure(face)/measure(ref. face),
3266  // det(J1)=measure(element)/measure(ref. element),
3267  // and the ratios measure(ref. element)/measure(ref. face) are
3268  // compatible for all element/face pairs.
3269  // For example: meas(ref. tetrahedron)/meas(ref. triangle) = 1/3, and
3270  // for any tetrahedron vol(tet)=(1/3)*height*area(base).
3271  // For interior faces: q_e/h_e=(q1/h1+q2/h2)/2.
3272 
3273  dshape1.Mult(nh, dshape1dn);
3274  for (int i = 0; i < ndof1; i++)
3275  for (int j = 0; j < ndof1; j++)
3276  {
3277  elmat(i, j) += shape1(i) * dshape1dn(j);
3278  }
3279 
3280  if (ndof2)
3281  {
3282  el2.CalcShape(eip2, shape2);
3283  el2.CalcDShape(eip2, dshape2);
3284  w = ip.weight/2/Trans.Elem2->Weight();
3285  if (!MQ)
3286  {
3287  if (Q)
3288  {
3289  w *= Q->Eval(*Trans.Elem2, eip2);
3290  }
3291  ni.Set(w, nor);
3292  }
3293  else
3294  {
3295  nh.Set(w, nor);
3296  MQ->Eval(mq, *Trans.Elem2, eip2);
3297  mq.MultTranspose(nh, ni);
3298  }
3299  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
3300  adjJ.Mult(ni, nh);
3301  if (kappa_is_nonzero)
3302  {
3303  wq += ni * nor;
3304  }
3305 
3306  dshape2.Mult(nh, dshape2dn);
3307 
3308  for (int i = 0; i < ndof1; i++)
3309  for (int j = 0; j < ndof2; j++)
3310  {
3311  elmat(i, ndof1 + j) += shape1(i) * dshape2dn(j);
3312  }
3313 
3314  for (int i = 0; i < ndof2; i++)
3315  for (int j = 0; j < ndof1; j++)
3316  {
3317  elmat(ndof1 + i, j) -= shape2(i) * dshape1dn(j);
3318  }
3319 
3320  for (int i = 0; i < ndof2; i++)
3321  for (int j = 0; j < ndof2; j++)
3322  {
3323  elmat(ndof1 + i, ndof1 + j) -= shape2(i) * dshape2dn(j);
3324  }
3325  }
3326 
3327  if (kappa_is_nonzero)
3328  {
3329  // only assemble the lower triangular part of jmat
3330  wq *= kappa;
3331  for (int i = 0; i < ndof1; i++)
3332  {
3333  const double wsi = wq*shape1(i);
3334  for (int j = 0; j <= i; j++)
3335  {
3336  jmat(i, j) += wsi * shape1(j);
3337  }
3338  }
3339  if (ndof2)
3340  {
3341  for (int i = 0; i < ndof2; i++)
3342  {
3343  const int i2 = ndof1 + i;
3344  const double wsi = wq*shape2(i);
3345  for (int j = 0; j < ndof1; j++)
3346  {
3347  jmat(i2, j) -= wsi * shape1(j);
3348  }
3349  for (int j = 0; j <= i; j++)
3350  {
3351  jmat(i2, ndof1 + j) += wsi * shape2(j);
3352  }
3353  }
3354  }
3355  }
3356  }
3357 
3358  // elmat := -elmat + sigma*elmat^t + jmat
3359  if (kappa_is_nonzero)
3360  {
3361  for (int i = 0; i < ndofs; i++)
3362  {
3363  for (int j = 0; j < i; j++)
3364  {
3365  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
3366  elmat(i,j) = sigma*aji - aij + mij;
3367  elmat(j,i) = sigma*aij - aji + mij;
3368  }
3369  elmat(i,i) = (sigma - 1.)*elmat(i,i) + jmat(i,i);
3370  }
3371  }
3372  else
3373  {
3374  for (int i = 0; i < ndofs; i++)
3375  {
3376  for (int j = 0; j < i; j++)
3377  {
3378  double aij = elmat(i,j), aji = elmat(j,i);
3379  elmat(i,j) = sigma*aji - aij;
3380  elmat(j,i) = sigma*aij - aji;
3381  }
3382  elmat(i,i) *= (sigma - 1.);
3383  }
3384  }
3385 }
3386 
3387 
3388 // static method
3389 void DGElasticityIntegrator::AssembleBlock(
3390  const int dim, const int row_ndofs, const int col_ndofs,
3391  const int row_offset, const int col_offset,
3392  const double jmatcoef, const Vector &col_nL, const Vector &col_nM,
3393  const Vector &row_shape, const Vector &col_shape,
3394  const Vector &col_dshape_dnM, const DenseMatrix &col_dshape,
3395  DenseMatrix &elmat, DenseMatrix &jmat)
3396 {
3397  for (int jm = 0, j = col_offset; jm < dim; ++jm)
3398  {
3399  for (int jdof = 0; jdof < col_ndofs; ++jdof, ++j)
3400  {
3401  const double t2 = col_dshape_dnM(jdof);
3402  for (int im = 0, i = row_offset; im < dim; ++im)
3403  {
3404  const double t1 = col_dshape(jdof, jm) * col_nL(im);
3405  const double t3 = col_dshape(jdof, im) * col_nM(jm);
3406  const double tt = t1 + ((im == jm) ? t2 : 0.0) + t3;
3407  for (int idof = 0; idof < row_ndofs; ++idof, ++i)
3408  {
3409  elmat(i, j) += row_shape(idof) * tt;
3410  }
3411  }
3412  }
3413  }
3414 
3415  if (jmatcoef == 0.0) { return; }
3416 
3417  for (int d = 0; d < dim; ++d)
3418  {
3419  const int jo = col_offset + d*col_ndofs;
3420  const int io = row_offset + d*row_ndofs;
3421  for (int jdof = 0, j = jo; jdof < col_ndofs; ++jdof, ++j)
3422  {
3423  const double sj = jmatcoef * col_shape(jdof);
3424  for (int i = max(io,j), idof = i - io; idof < row_ndofs; ++idof, ++i)
3425  {
3426  jmat(i, j) += row_shape(idof) * sj;
3427  }
3428  }
3429  }
3430 }
3431 
3432 void DGElasticityIntegrator::AssembleFaceMatrix(
3433  const FiniteElement &el1, const FiniteElement &el2,
3434  FaceElementTransformations &Trans, DenseMatrix &elmat)
3435 {
3436 #ifdef MFEM_THREAD_SAFE
3437  // For descriptions of these variables, see the class declaration.
3438  Vector shape1, shape2;
3439  DenseMatrix dshape1, dshape2;
3440  DenseMatrix adjJ;
3441  DenseMatrix dshape1_ps, dshape2_ps;
3442  Vector nor;
3443  Vector nL1, nL2;
3444  Vector nM1, nM2;
3445  Vector dshape1_dnM, dshape2_dnM;
3446  DenseMatrix jmat;
3447 #endif
3448 
3449  const int dim = el1.GetDim();
3450  const int ndofs1 = el1.GetDof();
3451  const int ndofs2 = (Trans.Elem2No >= 0) ? el2.GetDof() : 0;
3452  const int nvdofs = dim*(ndofs1 + ndofs2);
3453 
3454  // Initially 'elmat' corresponds to the term:
3455  // < { sigma(u) . n }, [v] > =
3456  // < { (lambda div(u) I + mu (grad(u) + grad(u)^T)) . n }, [v] >
3457  // But eventually, it's going to be replaced by:
3458  // elmat := -elmat + alpha*elmat^T + jmat
3459  elmat.SetSize(nvdofs);
3460  elmat = 0.;
3461 
3462  const bool kappa_is_nonzero = (kappa != 0.0);
3463  if (kappa_is_nonzero)
3464  {
3465  jmat.SetSize(nvdofs);
3466  jmat = 0.;
3467  }
3468 
3469  adjJ.SetSize(dim);
3470  shape1.SetSize(ndofs1);
3471  dshape1.SetSize(ndofs1, dim);
3472  dshape1_ps.SetSize(ndofs1, dim);
3473  nor.SetSize(dim);
3474  nL1.SetSize(dim);
3475  nM1.SetSize(dim);
3476  dshape1_dnM.SetSize(ndofs1);
3477 
3478  if (ndofs2)
3479  {
3480  shape2.SetSize(ndofs2);
3481  dshape2.SetSize(ndofs2, dim);
3482  dshape2_ps.SetSize(ndofs2, dim);
3483  nL2.SetSize(dim);
3484  nM2.SetSize(dim);
3485  dshape2_dnM.SetSize(ndofs2);
3486  }
3487 
3488  const IntegrationRule *ir = IntRule;
3489  if (ir == NULL)
3490  {
3491  // a simple choice for the integration order; is this OK?
3492  const int order = 2 * max(el1.GetOrder(), ndofs2 ? el2.GetOrder() : 0);
3493  ir = &IntRules.Get(Trans.GetGeometryType(), order);
3494  }
3495 
3496  for (int pind = 0; pind < ir->GetNPoints(); ++pind)
3497  {
3498  const IntegrationPoint &ip = ir->IntPoint(pind);
3499 
3500  // Set the integration point in the face and the neighboring elements
3501  Trans.SetAllIntPoints(&ip);
3502 
3503  // Access the neighboring elements' integration points
3504  // Note: eip2 will only contain valid data if Elem2 exists
3505  const IntegrationPoint &eip1 = Trans.GetElement1IntPoint();
3506  const IntegrationPoint &eip2 = Trans.GetElement2IntPoint();
3507 
3508  el1.CalcShape(eip1, shape1);
3509  el1.CalcDShape(eip1, dshape1);
3510 
3511  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
3512  Mult(dshape1, adjJ, dshape1_ps);
3513 
3514  if (dim == 1)
3515  {
3516  nor(0) = 2*eip1.x - 1.0;
3517  }
3518  else
3519  {
3520  CalcOrtho(Trans.Jacobian(), nor);
3521  }
3522 
3523  double w, wLM;
3524  if (ndofs2)
3525  {
3526  el2.CalcShape(eip2, shape2);
3527  el2.CalcDShape(eip2, dshape2);
3528  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
3529  Mult(dshape2, adjJ, dshape2_ps);
3530 
3531  w = ip.weight/2;
3532  const double w2 = w / Trans.Elem2->Weight();
3533  const double wL2 = w2 * lambda->Eval(*Trans.Elem2, eip2);
3534  const double wM2 = w2 * mu->Eval(*Trans.Elem2, eip2);
3535  nL2.Set(wL2, nor);
3536  nM2.Set(wM2, nor);
3537  wLM = (wL2 + 2.0*wM2);
3538  dshape2_ps.Mult(nM2, dshape2_dnM);
3539  }
3540  else
3541  {
3542  w = ip.weight;
3543  wLM = 0.0;
3544  }
3545 
3546  {
3547  const double w1 = w / Trans.Elem1->Weight();
3548  const double wL1 = w1 * lambda->Eval(*Trans.Elem1, eip1);
3549  const double wM1 = w1 * mu->Eval(*Trans.Elem1, eip1);
3550  nL1.Set(wL1, nor);
3551  nM1.Set(wM1, nor);
3552  wLM += (wL1 + 2.0*wM1);
3553  dshape1_ps.Mult(nM1, dshape1_dnM);
3554  }
3555 
3556  const double jmatcoef = kappa * (nor*nor) * wLM;
3557 
3558  // (1,1) block
3559  AssembleBlock(
3560  dim, ndofs1, ndofs1, 0, 0, jmatcoef, nL1, nM1,
3561  shape1, shape1, dshape1_dnM, dshape1_ps, elmat, jmat);
3562 
3563  if (ndofs2 == 0) { continue; }
3564 
3565  // In both elmat and jmat, shape2 appears only with a minus sign.
3566  shape2.Neg();
3567 
3568  // (1,2) block
3569  AssembleBlock(
3570  dim, ndofs1, ndofs2, 0, dim*ndofs1, jmatcoef, nL2, nM2,
3571  shape1, shape2, dshape2_dnM, dshape2_ps, elmat, jmat);
3572  // (2,1) block
3573  AssembleBlock(
3574  dim, ndofs2, ndofs1, dim*ndofs1, 0, jmatcoef, nL1, nM1,
3575  shape2, shape1, dshape1_dnM, dshape1_ps, elmat, jmat);
3576  // (2,2) block
3577  AssembleBlock(
3578  dim, ndofs2, ndofs2, dim*ndofs1, dim*ndofs1, jmatcoef, nL2, nM2,
3579  shape2, shape2, dshape2_dnM, dshape2_ps, elmat, jmat);
3580  }
3581 
3582  // elmat := -elmat + alpha*elmat^t + jmat
3583  if (kappa_is_nonzero)
3584  {
3585  for (int i = 0; i < nvdofs; ++i)
3586  {
3587  for (int j = 0; j < i; ++j)
3588  {
3589  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
3590  elmat(i,j) = alpha*aji - aij + mij;
3591  elmat(j,i) = alpha*aij - aji + mij;
3592  }
3593  elmat(i,i) = (alpha - 1.)*elmat(i,i) + jmat(i,i);
3594  }
3595  }
3596  else
3597  {
3598  for (int i = 0; i < nvdofs; ++i)
3599  {
3600  for (int j = 0; j < i; ++j)
3601  {
3602  double aij = elmat(i,j), aji = elmat(j,i);
3603  elmat(i,j) = alpha*aji - aij;
3604  elmat(j,i) = alpha*aij - aji;
3605  }
3606  elmat(i,i) *= (alpha - 1.);
3607  }
3608  }
3609 }
3610 
3611 
3612 void TraceJumpIntegrator::AssembleFaceMatrix(
3613  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
3614  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
3615  DenseMatrix &elmat)
3616 {
3617  int i, j, face_ndof, ndof1, ndof2;
3618  int order;
3619 
3620  double w;
3621 
3622  face_ndof = trial_face_fe.GetDof();
3623  ndof1 = test_fe1.GetDof();
3624 
3625  face_shape.SetSize(face_ndof);
3626  shape1.SetSize(ndof1);
3627 
3628  if (Trans.Elem2No >= 0)
3629  {
3630  ndof2 = test_fe2.GetDof();
3631  shape2.SetSize(ndof2);
3632  }
3633  else
3634  {
3635  ndof2 = 0;
3636  }
3637 
3638  elmat.SetSize(ndof1 + ndof2, face_ndof);
3639  elmat = 0.0;
3640 
3641  const IntegrationRule *ir = IntRule;
3642  if (ir == NULL)
3643  {
3644  if (Trans.Elem2No >= 0)
3645  {
3646  order = max(test_fe1.GetOrder(), test_fe2.GetOrder());
3647  }
3648  else
3649  {
3650  order = test_fe1.GetOrder();
3651  }
3652  order += trial_face_fe.GetOrder();
3653  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
3654  {
3655  order += Trans.OrderW();
3656  }
3657  ir = &IntRules.Get(Trans.GetGeometryType(), order);
3658  }
3659 
3660  for (int p = 0; p < ir->GetNPoints(); p++)
3661  {
3662  const IntegrationPoint &ip = ir->IntPoint(p);
3663 
3664  // Set the integration point in the face and the neighboring elements
3665  Trans.SetAllIntPoints(&ip);
3666 
3667  // Access the neighboring elements' integration points
3668  // Note: eip2 will only contain valid data if Elem2 exists
3669  const IntegrationPoint &eip1 = Trans.GetElement1IntPoint();
3670  const IntegrationPoint &eip2 = Trans.GetElement2IntPoint();
3671 
3672  // Trace finite element shape function
3673  trial_face_fe.CalcShape(ip, face_shape);
3674  // Side 1 finite element shape function
3675  test_fe1.CalcShape(eip1, shape1);
3676  if (ndof2)
3677  {
3678  // Side 2 finite element shape function
3679  test_fe2.CalcShape(eip2, shape2);
3680  }
3681  w = ip.weight;
3682  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
3683  {
3684  w *= Trans.Weight();
3685  }
3686  face_shape *= w;
3687  for (i = 0; i < ndof1; i++)
3688  for (j = 0; j < face_ndof; j++)
3689  {
3690  elmat(i, j) += shape1(i) * face_shape(j);
3691  }
3692  if (ndof2)
3693  {
3694  // Subtract contribution from side 2
3695  for (i = 0; i < ndof2; i++)
3696  for (j = 0; j < face_ndof; j++)
3697  {
3698  elmat(ndof1+i, j) -= shape2(i) * face_shape(j);
3699  }
3700  }
3701  }
3702 }
3703 
3704 void NormalTraceJumpIntegrator::AssembleFaceMatrix(
3705  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
3706  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
3707  DenseMatrix &elmat)
3708 {
3709  int i, j, face_ndof, ndof1, ndof2, dim;
3710  int order;
3711 
3712  MFEM_VERIFY(trial_face_fe.GetMapType() == FiniteElement::VALUE, "");
3713 
3714  face_ndof = trial_face_fe.GetDof();
3715  ndof1 = test_fe1.GetDof();
3716  dim = test_fe1.GetDim();
3717 
3718  face_shape.SetSize(face_ndof);
3719  normal.SetSize(dim);
3720  shape1.SetSize(ndof1,dim);
3721  shape1_n.SetSize(ndof1);
3722 
3723  if (Trans.Elem2No >= 0)
3724  {
3725  ndof2 = test_fe2.GetDof();
3726  shape2.SetSize(ndof2,dim);
3727  shape2_n.SetSize(ndof2);
3728  }
3729  else
3730  {
3731  ndof2 = 0;
3732  }
3733 
3734  elmat.SetSize(ndof1 + ndof2, face_ndof);
3735  elmat = 0.0;
3736 
3737  const IntegrationRule *ir = IntRule;
3738  if (ir == NULL)
3739  {
3740  if (Trans.Elem2No >= 0)
3741  {
3742  order = max(test_fe1.GetOrder(), test_fe2.GetOrder()) - 1;
3743  }
3744  else
3745  {
3746  order = test_fe1.GetOrder() - 1;
3747  }
3748  order += trial_face_fe.GetOrder();
3749  ir = &IntRules.Get(Trans.GetGeometryType(), order);
3750  }
3751 
3752  for (int p = 0; p < ir->GetNPoints(); p++)
3753  {
3754  const IntegrationPoint &ip = ir->IntPoint(p);
3755  IntegrationPoint eip1, eip2;
3756  // Trace finite element shape function
3757  trial_face_fe.CalcShape(ip, face_shape);
3758  Trans.Loc1.Transf.SetIntPoint(&ip);
3759  CalcOrtho(Trans.Loc1.Transf.Jacobian(), normal);
3760  // Side 1 finite element shape function
3761  Trans.Loc1.Transform(ip, eip1);
3762  test_fe1.CalcVShape(eip1, shape1);
3763  shape1.Mult(normal, shape1_n);
3764  if (ndof2)
3765  {
3766  // Side 2 finite element shape function
3767  Trans.Loc2.Transform(ip, eip2);
3768  test_fe2.CalcVShape(eip2, shape2);
3769  Trans.Loc2.Transf.SetIntPoint(&ip);
3770  CalcOrtho(Trans.Loc2.Transf.Jacobian(), normal);
3771  shape2.Mult(normal, shape2_n);
3772  }
3773  face_shape *= ip.weight;
3774  for (i = 0; i < ndof1; i++)
3775  for (j = 0; j < face_ndof; j++)
3776  {
3777  elmat(i, j) -= shape1_n(i) * face_shape(j);
3778  }
3779  if (ndof2)
3780  {
3781  // Subtract contribution from side 2
3782  for (i = 0; i < ndof2; i++)
3783  for (j = 0; j < face_ndof; j++)
3784  {
3785  elmat(ndof1+i, j) += shape2_n(i) * face_shape(j);
3786  }
3787  }
3788  }
3789 }
3790 
3791 
3792 void NormalInterpolator::AssembleElementMatrix2(
3793  const FiniteElement &dom_fe, const FiniteElement &ran_fe,
3794  ElementTransformation &Trans, DenseMatrix &elmat)
3795 {
3796  int spaceDim = Trans.GetSpaceDim();
3797  elmat.SetSize(ran_fe.GetDof(), spaceDim*dom_fe.GetDof());
3798  Vector n(spaceDim), shape(dom_fe.GetDof());
3799 
3800  const IntegrationRule &ran_nodes = ran_fe.GetNodes();
3801  for (int i = 0; i < ran_nodes.Size(); i++)
3802  {
3803  const IntegrationPoint &ip = ran_nodes.IntPoint(i);
3804  Trans.SetIntPoint(&ip);
3805  CalcOrtho(Trans.Jacobian(), n);
3806  dom_fe.CalcShape(ip, shape);
3807  for (int j = 0; j < shape.Size(); j++)
3808  {
3809  for (int d = 0; d < spaceDim; d++)
3810  {
3811  elmat(i, j+d*shape.Size()) = shape(j)*n(d);
3812  }
3813  }
3814  }
3815 }
3816 
3817 
3818 namespace internal
3819 {
3820 
3821 // Scalar shape functions scaled by scalar coefficient.
3822 // Used in the implementation of class ScalarProductInterpolator below.
3823 struct ShapeCoefficient : public VectorCoefficient
3824 {
3825  Coefficient &Q;
3826  const FiniteElement &fe;
3827 
3828  ShapeCoefficient(Coefficient &q, const FiniteElement &fe_)
3829  : VectorCoefficient(fe_.GetDof()), Q(q), fe(fe_) { }
3830 
3831  using VectorCoefficient::Eval;
3832  virtual void Eval(Vector &V, ElementTransformation &T,
3833  const IntegrationPoint &ip)
3834  {
3835  V.SetSize(vdim);
3836  fe.CalcPhysShape(T, V);
3837  V *= Q.Eval(T, ip);
3838  }
3839 };
3840 
3841 }
3842 
3843 void
3844 ScalarProductInterpolator::AssembleElementMatrix2(const FiniteElement &dom_fe,
3845  const FiniteElement &ran_fe,
3846  ElementTransformation &Trans,
3847  DenseMatrix &elmat)
3848 {
3849  internal::ShapeCoefficient dom_shape_coeff(*Q, dom_fe);
3850 
3851  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
3852 
3853  Vector elmat_as_vec(elmat.Data(), ran_fe.GetDof()*dom_fe.GetDof());
3854 
3855  ran_fe.Project(dom_shape_coeff, Trans, elmat_as_vec);
3856 }
3857 
3858 
3859 void
3860 ScalarVectorProductInterpolator::AssembleElementMatrix2(
3861  const FiniteElement &dom_fe,
3862  const FiniteElement &ran_fe,
3863  ElementTransformation &Trans,
3864  DenseMatrix &elmat)
3865 {
3866  // Vector shape functions scaled by scalar coefficient
3867  struct VShapeCoefficient : public MatrixCoefficient
3868  {
3869  Coefficient &Q;
3870  const FiniteElement &fe;
3871 
3872  VShapeCoefficient(Coefficient &q, const FiniteElement &fe_, int sdim)
3873  : MatrixCoefficient(fe_.GetDof(), sdim), Q(q), fe(fe_) { }
3874 
3875  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
3876  const IntegrationPoint &ip)
3877  {
3878  M.SetSize(height, width);
3879  fe.CalcPhysVShape(T, M);
3880  M *= Q.Eval(T, ip);
3881  }
3882  };
3883 
3884  VShapeCoefficient dom_shape_coeff(*Q, dom_fe, Trans.GetSpaceDim());
3885 
3886  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
3887 
3888  Vector elmat_as_vec(elmat.Data(), ran_fe.GetDof()*dom_fe.GetDof());
3889 
3890  ran_fe.ProjectMatrixCoefficient(dom_shape_coeff, Trans, elmat_as_vec);
3891 }
3892 
3893 
3894 void
3895 VectorScalarProductInterpolator::AssembleElementMatrix2(
3896  const FiniteElement &dom_fe,
3897  const FiniteElement &ran_fe,
3898  ElementTransformation &Trans,
3899  DenseMatrix &elmat)
3900 {
3901  // Scalar shape functions scaled by vector coefficient
3902  struct VecShapeCoefficient : public MatrixCoefficient
3903  {
3904  VectorCoefficient &VQ;
3905  const FiniteElement &fe;
3906  Vector vc, shape;
3907 
3908  VecShapeCoefficient(VectorCoefficient &vq, const FiniteElement &fe_)
3909  : MatrixCoefficient(fe_.GetDof(), vq.GetVDim()), VQ(vq), fe(fe_),
3910  vc(width), shape(height) { }
3911 
3912  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
3913  const IntegrationPoint &ip)
3914  {
3915  M.SetSize(height, width);
3916  VQ.Eval(vc, T, ip);
3917  fe.CalcPhysShape(T, shape);
3918  MultVWt(shape, vc, M);
3919  }
3920  };
3921 
3922  VecShapeCoefficient dom_shape_coeff(*VQ, dom_fe);
3923 
3924  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
3925 
3926  Vector elmat_as_vec(elmat.Data(), ran_fe.GetDof()*dom_fe.GetDof());
3927 
3928  ran_fe.ProjectMatrixCoefficient(dom_shape_coeff, Trans, elmat_as_vec);
3929 }
3930 
3931 
3932 void
3933 ScalarCrossProductInterpolator::AssembleElementMatrix2(
3934  const FiniteElement &dom_fe,
3935  const FiniteElement &ran_fe,
3936  ElementTransformation &Trans,
3937  DenseMatrix &elmat)
3938 {
3939  // Vector coefficient product with vector shape functions
3940  struct VCrossVShapeCoefficient : public VectorCoefficient
3941  {
3942  VectorCoefficient &VQ;
3943  const FiniteElement &fe;
3944  DenseMatrix vshape;
3945  Vector vc;
3946 
3947  VCrossVShapeCoefficient(VectorCoefficient &vq, const FiniteElement &fe_)
3948  : VectorCoefficient(fe_.GetDof()), VQ(vq), fe(fe_),
3949  vshape(vdim, vq.GetVDim()), vc(vq.GetVDim()) { }
3950 
3951  using VectorCoefficient::Eval;
3952  virtual void Eval(Vector &V, ElementTransformation &T,
3953  const IntegrationPoint &ip)
3954  {
3955  V.SetSize(vdim);
3956  VQ.Eval(vc, T, ip);
3957  fe.CalcPhysVShape(T, vshape);
3958  for (int k = 0; k < vdim; k++)
3959  {
3960  V(k) = vc(0) * vshape(k,1) - vc(1) * vshape(k,0);
3961  }
3962  }
3963  };
3964 
3965  VCrossVShapeCoefficient dom_shape_coeff(*VQ, dom_fe);
3966 
3967  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
3968 
3969  Vector elmat_as_vec(elmat.Data(), elmat.Height()*elmat.Width());
3970 
3971  ran_fe.Project(dom_shape_coeff, Trans, elmat_as_vec);
3972 }
3973 
3974 void
3975 VectorCrossProductInterpolator::AssembleElementMatrix2(
3976  const FiniteElement &dom_fe,
3977  const FiniteElement &ran_fe,
3978  ElementTransformation &Trans,
3979  DenseMatrix &elmat)
3980 {
3981  // Vector coefficient product with vector shape functions
3982  struct VCrossVShapeCoefficient : public MatrixCoefficient
3983  {
3984  VectorCoefficient &VQ;
3985  const FiniteElement &fe;
3986  DenseMatrix vshape;
3987  Vector vc;
3988 
3989  VCrossVShapeCoefficient(VectorCoefficient &vq, const FiniteElement &fe_)
3990  : MatrixCoefficient(fe_.GetDof(), vq.GetVDim()), VQ(vq), fe(fe_),
3991  vshape(height, width), vc(width)
3992  {
3993  MFEM_ASSERT(width == 3, "");
3994  }
3995 
3996  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
3997  const IntegrationPoint &ip)
3998  {
3999  M.SetSize(height, width);
4000  VQ.Eval(vc, T, ip);
4001  fe.CalcPhysVShape(T, vshape);
4002  for (int k = 0; k < height; k++)
4003  {
4004  M(k,0) = vc(1) * vshape(k,2) - vc(2) * vshape(k,1);
4005  M(k,1) = vc(2) * vshape(k,0) - vc(0) * vshape(k,2);
4006  M(k,2) = vc(0) * vshape(k,1) - vc(1) * vshape(k,0);
4007  }
4008  }
4009  };
4010 
4011  VCrossVShapeCoefficient dom_shape_coeff(*VQ, dom_fe);
4012 
4013  if (ran_fe.GetRangeType() == FiniteElement::SCALAR)
4014  {
4015  elmat.SetSize(ran_fe.GetDof()*VQ->GetVDim(),dom_fe.GetDof());
4016  }
4017  else
4018  {
4019  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
4020  }
4021 
4022  Vector elmat_as_vec(elmat.Data(), elmat.Height()*elmat.Width());
4023 
4024  ran_fe.ProjectMatrixCoefficient(dom_shape_coeff, Trans, elmat_as_vec);
4025 }
4026 
4027 
4028 namespace internal
4029 {
4030 
4031 // Vector shape functions dot product with a vector coefficient.
4032 // Used in the implementation of class VectorInnerProductInterpolator below.
4033 struct VDotVShapeCoefficient : public VectorCoefficient
4034 {
4035  VectorCoefficient &VQ;
4036  const FiniteElement &fe;
4037  DenseMatrix vshape;
4038  Vector vc;
4039 
4040  VDotVShapeCoefficient(VectorCoefficient &vq, const FiniteElement &fe_)
4041  : VectorCoefficient(fe_.GetDof()), VQ(vq), fe(fe_),
4042  vshape(vdim, vq.GetVDim()), vc(vq.GetVDim()) { }
4043 
4044  using VectorCoefficient::Eval;
4045  virtual void Eval(Vector &V, ElementTransformation &T,
4046  const IntegrationPoint &ip)
4047  {
4048  V.SetSize(vdim);
4049  VQ.Eval(vc, T, ip);
4050  fe.CalcPhysVShape(T, vshape);
4051  vshape.Mult(vc, V);
4052  }
4053 };
4054 
4055 }
4056 
4057 void
4058 VectorInnerProductInterpolator::AssembleElementMatrix2(
4059  const FiniteElement &dom_fe,
4060  const FiniteElement &ran_fe,
4061  ElementTransformation &Trans,
4062  DenseMatrix &elmat)
4063 {
4064  internal::VDotVShapeCoefficient dom_shape_coeff(*VQ, dom_fe);
4065 
4066  elmat.SetSize(ran_fe.GetDof(),dom_fe.GetDof());
4067 
4068  Vector elmat_as_vec(elmat.Data(), elmat.Height()*elmat.Width());
4069 
4070  ran_fe.Project(dom_shape_coeff, Trans, elmat_as_vec);
4071 }
4072 
4073 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
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:2376
const DenseMatrix & AdjugateJacobian()
Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoin...
Definition: eltrans.hpp:135
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:2761
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual int OrderGrad(const FiniteElement *fe) const =0
Return the order of .
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:162
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_base.cpp:39
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:2742
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
Base class for vector Coefficients that optionally depend on time and space.
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition: fe_base.cpp:125
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:284
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2074
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:524
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
void AddMult_a_ABt(double a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
Definition: densemat.cpp:2574
const DenseMatrix & InverseJacobian()
Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:145
double kappa
Definition: ex24.cpp:54
const IntegrationPoint & GetElement2IntPoint()
Get a const reference to the integration point in neighboring element 2 corresponding to the currentl...
Definition: eltrans.hpp:587
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe_base.hpp:349
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:115
ElementTransformation * Elem2
Definition: eltrans.hpp:522
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2288
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
int GetCurlDim() const
Returns the dimension of the curl for vector-valued finite elements.
Definition: fe_base.hpp:317
void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
Definition: densemat.cpp:2715
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:300
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
Definition: fe_base.cpp:70
virtual int OrderW() const
Return the order of the determinant of the Jacobian (weight) of the transformation.
Definition: eltrans.cpp:457
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:167
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:188
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:90
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:524
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:2806
double b
Definition: lissajous.cpp:42
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition: eltrans.hpp:566
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
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...
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:618
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:2828
int GetVDim()
Returns dimension of the vector.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2182
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:1359
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:2731
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:2481
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
virtual int Order() const =0
Return the order of the current element we are using for the transformation.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition: fe_base.cpp:181
void ClearExternalData()
Definition: densemat.hpp:95
Base class for Matrix Coefficients that optionally depend on time and space.
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:546
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:2314
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:285
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_base.cpp:51
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:1648
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition: eltrans.hpp:577
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:156
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:267
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
double a
Definition: lissajous.cpp:41
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_base.cpp:64
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:2538
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation (&quot;projection&quot;) in the local...
Definition: fe_base.cpp:143
double mu
Definition: ex25.cpp:139
int dim
Definition: ex24.cpp:53
IsoparametricTransformation Transf
Definition: eltrans.hpp:464
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:1417
ElementTransformation * Elem1
Definition: eltrans.hpp:522
const double alpha
Definition: ex15.cpp:369
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:2633
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:160
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:382
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
Definition: densemat.cpp:2332
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
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
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition: fe_base.hpp:340
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:2690
void Neg()
(*this) = -(*this)
Definition: vector.cpp:292
double sigma(const Vector &x)
Definition: maxwell.cpp:102