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