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