MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Implementation of Bilinear Form Integrators
13 
14 #include "fem.hpp"
15 #include <cmath>
16 #include <algorithm>
17 
18 using namespace std;
19 
20 namespace mfem
21 {
22 
23 void BilinearFormIntegrator::AssembleElementMatrix (
25  DenseMatrix &elmat )
26 {
27  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix (...)\n"
28  " is not implemented fot this class.");
29 }
30 
31 void BilinearFormIntegrator::AssembleElementMatrix2 (
32  const FiniteElement &el1, const FiniteElement &el2,
34 {
35  mfem_error ("BilinearFormIntegrator::AssembleElementMatrix2 (...)\n"
36  " is not implemented fot this class.");
37 }
38 
39 void BilinearFormIntegrator::AssembleFaceMatrix (
40  const FiniteElement &el1, const FiniteElement &el2,
42 {
43  mfem_error ("BilinearFormIntegrator::AssembleFaceMatrix (...)\n"
44  " is not implemented fot this class.");
45 }
46 
47 void BilinearFormIntegrator::AssembleFaceMatrix(
48  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
50  DenseMatrix &elmat)
51 {
52  MFEM_ABORT("AssembleFaceMatrix (mixed form) is not implemented for this"
53  " Integrator class.");
54 }
55 
56 void BilinearFormIntegrator::AssembleElementVector(
57  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
58  Vector &elvect)
59 {
60  mfem_error("BilinearFormIntegrator::AssembleElementVector\n"
61  " is not implemented fot this class.");
62 }
63 
64 
65 void TransposeIntegrator::AssembleElementMatrix (
67 {
68  bfi -> AssembleElementMatrix (el, Trans, bfi_elmat);
69  // elmat = bfi_elmat^t
70  elmat.Transpose (bfi_elmat);
71 }
72 
73 void TransposeIntegrator::AssembleElementMatrix2 (
74  const FiniteElement &trial_fe, const FiniteElement &test_fe,
76 {
77  bfi -> AssembleElementMatrix2 (test_fe, trial_fe, Trans, bfi_elmat);
78  // elmat = bfi_elmat^t
79  elmat.Transpose (bfi_elmat);
80 }
81 
82 void TransposeIntegrator::AssembleFaceMatrix (
83  const FiniteElement &el1, const FiniteElement &el2,
85 {
86  bfi -> AssembleFaceMatrix (el1, el2, Trans, bfi_elmat);
87  // elmat = bfi_elmat^t
88  elmat.Transpose (bfi_elmat);
89 }
90 
91 void LumpedIntegrator::AssembleElementMatrix (
93 {
94  bfi -> AssembleElementMatrix (el, Trans, elmat);
95  elmat.Lump();
96 }
97 
98 void InverseIntegrator::AssembleElementMatrix(
100 {
101  integrator->AssembleElementMatrix(el, Trans, elmat);
102  elmat.Invert();
103 }
104 
105 void SumIntegrator::AssembleElementMatrix(
107 {
108  MFEM_ASSERT(integrators.Size() > 0, "empty SumIntegrator.");
109 
110  integrators[0]->AssembleElementMatrix(el, Trans, elmat);
111  for (int i = 1; i < integrators.Size(); i++)
112  {
113  integrators[i]->AssembleElementMatrix(el, Trans, elem_mat);
114  elmat += elem_mat;
115  }
116 }
117 
118 SumIntegrator::~SumIntegrator()
119 {
120  if (own_integrators)
121  {
122  for (int i = 0; i < integrators.Size(); i++)
123  {
124  delete integrators[i];
125  }
126  }
127 }
128 
129 
130 void DiffusionIntegrator::AssembleElementMatrix
132  DenseMatrix &elmat )
133 {
134  int nd = el.GetDof();
135  int dim = el.GetDim();
136  int spaceDim = Trans.GetSpaceDim();
137  bool square = (dim == spaceDim);
138  double w;
139 
140 #ifdef MFEM_THREAD_SAFE
141  DenseMatrix dshape(nd,dim), dshapedxt(nd,spaceDim), invdfdx(dim,spaceDim);
142 #else
143  dshape.SetSize(nd,dim);
144  dshapedxt.SetSize(nd,spaceDim);
145  invdfdx.SetSize(dim,spaceDim);
146 #endif
147  elmat.SetSize(nd);
148 
149  const IntegrationRule *ir = IntRule;
150  if (ir == NULL)
151  {
152  int order;
153  if (el.Space() == FunctionSpace::Pk)
154  {
155  order = 2*el.GetOrder() - 2;
156  }
157  else
158  // order = 2*el.GetOrder() - 2; // <-- this seems to work fine too
159  {
160  order = 2*el.GetOrder() + dim - 1;
161  }
162 
163  if (el.Space() == FunctionSpace::rQk)
164  {
165  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
166  }
167  else
168  {
169  ir = &IntRules.Get(el.GetGeomType(), order);
170  }
171  }
172 
173  elmat = 0.0;
174  for (int i = 0; i < ir->GetNPoints(); i++)
175  {
176  const IntegrationPoint &ip = ir->IntPoint(i);
177  el.CalcDShape(ip, dshape);
178 
179  Trans.SetIntPoint(&ip);
180  // Compute invdfdx = / adj(J), if J is square
181  // \ adj(J^t.J).J^t, otherwise
182  CalcAdjugate(Trans.Jacobian(), invdfdx);
183  w = Trans.Weight();
184  w = ip.weight / (square ? w : w*w*w);
185  Mult(dshape, invdfdx, dshapedxt);
186  if (!MQ)
187  {
188  if (Q)
189  {
190  w *= Q->Eval(Trans, ip);
191  }
192  AddMult_a_AAt(w, dshapedxt, elmat);
193  }
194  else
195  {
196  MQ->Eval(invdfdx, Trans, ip);
197  invdfdx *= w;
198  Mult(dshapedxt, invdfdx, dshape);
199  AddMultABt(dshape, dshapedxt, elmat);
200  }
201  }
202 }
203 
204 void DiffusionIntegrator::AssembleElementMatrix2(
205  const FiniteElement &trial_fe, const FiniteElement &test_fe,
206  ElementTransformation &Trans, DenseMatrix &elmat)
207 {
208  int tr_nd = trial_fe.GetDof();
209  int te_nd = test_fe.GetDof();
210  int dim = trial_fe.GetDim();
211  int spaceDim = Trans.GetSpaceDim();
212  bool square = (dim == spaceDim);
213  double w;
214 
215 #ifdef MFEM_THREAD_SAFE
216  DenseMatrix dshape(tr_nd, dim), dshapedxt(tr_nd, spaceDim);
217  DenseMatrix te_dshape(te_nd, dim), te_dshapedxt(te_nd, spaceDim);
218  DenseMatrix invdfdx(dim, spaceDim);
219 #else
220  dshape.SetSize(tr_nd, dim);
221  dshapedxt.SetSize(tr_nd, spaceDim);
222  te_dshape.SetSize(te_nd, dim);
223  te_dshapedxt.SetSize(te_nd, spaceDim);
224  invdfdx.SetSize(dim, spaceDim);
225 #endif
226  elmat.SetSize(te_nd, tr_nd);
227 
228  const IntegrationRule *ir = IntRule;
229  if (ir == NULL)
230  {
231  int order;
232  if (trial_fe.Space() == FunctionSpace::Pk)
233  {
234  order = trial_fe.GetOrder() + test_fe.GetOrder() - 2;
235  }
236  else
237  {
238  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim - 1;
239  }
240 
241  if (trial_fe.Space() == FunctionSpace::rQk)
242  {
243  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
244  }
245  else
246  {
247  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
248  }
249  }
250 
251  elmat = 0.0;
252  for (int i = 0; i < ir->GetNPoints(); i++)
253  {
254  const IntegrationPoint &ip = ir->IntPoint(i);
255  trial_fe.CalcDShape(ip, dshape);
256  test_fe.CalcDShape(ip, te_dshape);
257 
258  Trans.SetIntPoint(&ip);
259  CalcAdjugate(Trans.Jacobian(), invdfdx);
260  w = Trans.Weight();
261  w = ip.weight / (square ? w : w*w*w);
262  Mult(dshape, invdfdx, dshapedxt);
263  Mult(te_dshape, invdfdx, te_dshapedxt);
264  // invdfdx, dshape, and te_dshape no longer needed
265  if (!MQ)
266  {
267  if (Q)
268  {
269  w *= Q->Eval(Trans, ip);
270  }
271  dshapedxt *= w;
272  AddMultABt(te_dshapedxt, dshapedxt, elmat);
273  }
274  else
275  {
276  MQ->Eval(invdfdx, Trans, ip);
277  invdfdx *= w;
278  Mult(te_dshapedxt, invdfdx, te_dshape);
279  AddMultABt(te_dshape, dshapedxt, elmat);
280  }
281  }
282 }
283 
284 void DiffusionIntegrator::AssembleElementVector(
285  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
286  Vector &elvect)
287 {
288  int nd = el.GetDof();
289  int dim = el.GetDim();
290  double w;
291 
292 #ifdef MFEM_THREAD_SAFE
293  DenseMatrix dshape(nd,dim), invdfdx(dim), mq(dim);
294 #else
295  dshape.SetSize(nd,dim);
296  invdfdx.SetSize(dim);
297  mq.SetSize(dim);
298 #endif
299  vec.SetSize(dim);
300  pointflux.SetSize(dim);
301 
302  elvect.SetSize(nd);
303 
304  const IntegrationRule *ir = IntRule;
305  if (ir == NULL)
306  {
307  int order;
308  if (el.Space() == FunctionSpace::Pk)
309  {
310  order = 2*el.GetOrder() - 2;
311  }
312  else
313  // order = 2*el.GetOrder() - 2; // <-- this seems to work fine too
314  {
315  order = 2*el.GetOrder() + dim - 1;
316  }
317 
318  if (el.Space() == FunctionSpace::rQk)
319  {
320  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
321  }
322  else
323  {
324  ir = &IntRules.Get(el.GetGeomType(), order);
325  }
326  }
327 
328  elvect = 0.0;
329  for (int i = 0; i < ir->GetNPoints(); i++)
330  {
331  const IntegrationPoint &ip = ir->IntPoint(i);
332  el.CalcDShape(ip, dshape);
333 
334  Tr.SetIntPoint(&ip);
335  CalcAdjugate(Tr.Jacobian(), invdfdx); // invdfdx = adj(J)
336  w = ip.weight / Tr.Weight();
337 
338  if (!MQ)
339  {
340  dshape.MultTranspose(elfun, vec);
341  invdfdx.MultTranspose(vec, pointflux);
342  if (Q)
343  {
344  w *= Q->Eval(Tr, ip);
345  }
346  }
347  else
348  {
349 
350  dshape.MultTranspose(elfun, pointflux);
351  invdfdx.MultTranspose(pointflux, vec);
352  MQ->Eval(mq, Tr, ip);
353  mq.Mult(vec, pointflux);
354  }
355  pointflux *= w;
356  invdfdx.Mult(pointflux, vec);
357  dshape.AddMult(vec, elvect);
358  }
359 }
360 
361 void DiffusionIntegrator::ComputeElementFlux
363  Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef )
364 {
365  int i, j, nd, dim, spaceDim, fnd;
366 
367  nd = el.GetDof();
368  dim = el.GetDim();
369  spaceDim = Trans.GetSpaceDim();
370 
371 #ifdef MFEM_THREAD_SAFE
372  DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim);
373 #else
374  dshape.SetSize(nd,dim);
375  invdfdx.SetSize(dim, spaceDim);
376 #endif
377  vec.SetSize(dim);
378  pointflux.SetSize(spaceDim);
379 
380  const IntegrationRule &ir = fluxelem.GetNodes();
381  fnd = ir.GetNPoints();
382  flux.SetSize( fnd * spaceDim );
383 
384  for (i = 0; i < fnd; i++)
385  {
386  const IntegrationPoint &ip = ir.IntPoint(i);
387  el.CalcDShape(ip, dshape);
388  dshape.MultTranspose(u, vec);
389 
390  Trans.SetIntPoint (&ip);
391  CalcInverse(Trans.Jacobian(), invdfdx);
392  invdfdx.MultTranspose(vec, pointflux);
393 
394  if (!MQ)
395  {
396  if (Q && with_coef)
397  {
398  pointflux *= Q->Eval(Trans,ip);
399  }
400  for (j = 0; j < spaceDim; j++)
401  {
402  flux(fnd*j+i) = pointflux(j);
403  }
404  }
405  else
406  {
407  // assuming dim == spaceDim
408  MFEM_ASSERT(dim == spaceDim, "TODO");
409  MQ->Eval(invdfdx, Trans, ip);
410  invdfdx.Mult(pointflux, vec);
411  for (j = 0; j < dim; j++)
412  {
413  flux(fnd*j+i) = vec(j);
414  }
415  }
416  }
417 }
418 
419 double DiffusionIntegrator::ComputeFluxEnergy
420 ( const FiniteElement &fluxelem, ElementTransformation &Trans,
421  Vector &flux, Vector* d_energy)
422 {
423  int nd = fluxelem.GetDof();
424  int dim = fluxelem.GetDim();
425  int spaceDim = Trans.GetSpaceDim();
426 
427 #ifdef MFEM_THREAD_SAFE
428  DenseMatrix mq;
429 #endif
430 
431  shape.SetSize(nd);
432  pointflux.SetSize(spaceDim);
433  if (d_energy) { vec.SetSize(dim); }
434  if (MQ) { mq.SetSize(dim); }
435 
436  int order = 2 * fluxelem.GetOrder(); // <--
437  const IntegrationRule *ir = &IntRules.Get(fluxelem.GetGeomType(), order);
438 
439  double energy = 0.0;
440  if (d_energy) { *d_energy = 0.0; }
441 
442  for (int i = 0; i < ir->GetNPoints(); i++)
443  {
444  const IntegrationPoint &ip = ir->IntPoint(i);
445  fluxelem.CalcShape(ip, shape);
446 
447  pointflux = 0.0;
448  for (int k = 0; k < spaceDim; k++)
449  {
450  for (int j = 0; j < nd; j++)
451  {
452  pointflux(k) += flux(k*nd+j)*shape(j);
453  }
454  }
455 
456  Trans.SetIntPoint(&ip);
457  double w = Trans.Weight() * ip.weight;
458 
459  if (!MQ)
460  {
461  double e = (pointflux * pointflux);
462  if (Q) { e *= Q->Eval(Trans, ip); }
463  energy += w * e;
464  }
465  else
466  {
467  MQ->Eval(mq, Trans, ip);
468  energy += w * mq.InnerProduct(pointflux, pointflux);
469  }
470 
471  if (d_energy)
472  {
473  // transform pointflux to the ref. domain and integrate the components
474  Trans.Jacobian().MultTranspose(pointflux, vec);
475  for (int k = 0; k < dim; k++)
476  {
477  (*d_energy)[k] += w * vec[k] * vec[k];
478  }
479  // TODO: Q, MQ
480  }
481  }
482 
483  return energy;
484 }
485 
486 
487 void MassIntegrator::AssembleElementMatrix
489  DenseMatrix &elmat )
490 {
491  int nd = el.GetDof();
492  // int dim = el.GetDim();
493  double w;
494 
495  elmat.SetSize(nd);
496  shape.SetSize(nd);
497 
498  const IntegrationRule *ir = IntRule;
499  if (ir == NULL)
500  {
501  // int order = 2 * el.GetOrder();
502  int order = 2 * el.GetOrder() + Trans.OrderW();
503 
504  if (el.Space() == FunctionSpace::rQk)
505  {
506  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
507  }
508  else
509  {
510  ir = &IntRules.Get(el.GetGeomType(), order);
511  }
512  }
513 
514  elmat = 0.0;
515  for (int i = 0; i < ir->GetNPoints(); i++)
516  {
517  const IntegrationPoint &ip = ir->IntPoint(i);
518  el.CalcShape(ip, shape);
519 
520  Trans.SetIntPoint (&ip);
521  w = Trans.Weight() * ip.weight;
522  if (Q)
523  {
524  w *= Q -> Eval(Trans, ip);
525  }
526 
527  AddMult_a_VVt(w, shape, elmat);
528  }
529 }
530 
531 void MassIntegrator::AssembleElementMatrix2(
532  const FiniteElement &trial_fe, const FiniteElement &test_fe,
533  ElementTransformation &Trans, DenseMatrix &elmat)
534 {
535  int tr_nd = trial_fe.GetDof();
536  int te_nd = test_fe.GetDof();
537  // int dim = trial_fe.GetDim();
538  double w;
539 
540  elmat.SetSize (te_nd, tr_nd);
541  shape.SetSize (tr_nd);
542  te_shape.SetSize (te_nd);
543 
544  const IntegrationRule *ir = IntRule;
545  if (ir == NULL)
546  {
547  int order = trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW();
548 
549  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
550  }
551 
552  elmat = 0.0;
553  for (int i = 0; i < ir->GetNPoints(); i++)
554  {
555  const IntegrationPoint &ip = ir->IntPoint(i);
556  trial_fe.CalcShape(ip, shape);
557  test_fe.CalcShape(ip, te_shape);
558 
559  Trans.SetIntPoint (&ip);
560  w = Trans.Weight() * ip.weight;
561  if (Q)
562  {
563  w *= Q -> Eval(Trans, ip);
564  }
565 
566  te_shape *= w;
567  AddMultVWt(te_shape, shape, elmat);
568  }
569 }
570 
571 
572 void BoundaryMassIntegrator::AssembleFaceMatrix(
573  const FiniteElement &el1, const FiniteElement &el2,
575 {
576  MFEM_ASSERT(Trans.Elem2No < 0,
577  "support for interior faces is not implemented");
578 
579  int nd1 = el1.GetDof();
580  double w;
581 
582  elmat.SetSize(nd1);
583  shape.SetSize(nd1);
584 
585  const IntegrationRule *ir = IntRule;
586  if (ir == NULL)
587  {
588  int order = 2 * el1.GetOrder();
589 
590  ir = &IntRules.Get(Trans.FaceGeom, order);
591  }
592 
593  elmat = 0.0;
594  for (int i = 0; i < ir->GetNPoints(); i++)
595  {
596  const IntegrationPoint &ip = ir->IntPoint(i);
597  IntegrationPoint eip;
598  Trans.Loc1.Transform(ip, eip);
599  el1.CalcShape(eip, shape);
600 
601  Trans.Face->SetIntPoint(&ip);
602  w = Trans.Face->Weight() * ip.weight;
603  if (Q)
604  {
605  w *= Q -> Eval(*Trans.Face, ip);
606  }
607 
608  AddMult_a_VVt(w, shape, elmat);
609  }
610 }
611 
612 
613 void ConvectionIntegrator::AssembleElementMatrix(
614  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
615 {
616  int nd = el.GetDof();
617  int dim = el.GetDim();
618 
619  elmat.SetSize(nd);
620  dshape.SetSize(nd,dim);
621  adjJ.SetSize(dim);
622  shape.SetSize(nd);
623  vec2.SetSize(dim);
624  BdFidxT.SetSize(nd);
625 
626  Vector vec1;
627 
628  const IntegrationRule *ir = IntRule;
629  if (ir == NULL)
630  {
631  int order = Trans.OrderGrad(&el) + Trans.Order() + el.GetOrder();
632  ir = &IntRules.Get(el.GetGeomType(), order);
633  }
634 
635  Q.Eval(Q_ir, Trans, *ir);
636 
637  elmat = 0.0;
638  for (int i = 0; i < ir->GetNPoints(); i++)
639  {
640  const IntegrationPoint &ip = ir->IntPoint(i);
641  el.CalcDShape(ip, dshape);
642  el.CalcShape(ip, shape);
643 
644  Trans.SetIntPoint(&ip);
645  CalcAdjugate(Trans.Jacobian(), adjJ);
646  Q_ir.GetColumnReference(i, vec1);
647  vec1 *= alpha * ip.weight;
648 
649  adjJ.Mult(vec1, vec2);
650  dshape.Mult(vec2, BdFidxT);
651 
652  AddMultVWt(shape, BdFidxT, elmat);
653  }
654 }
655 
656 
657 void GroupConvectionIntegrator::AssembleElementMatrix(
658  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
659 {
660  int nd = el.GetDof();
661  int dim = el.GetDim();
662 
663  elmat.SetSize(nd);
664  dshape.SetSize(nd,dim);
665  adjJ.SetSize(dim);
666  shape.SetSize(nd);
667  grad.SetSize(nd,dim);
668 
669  const IntegrationRule *ir = IntRule;
670  if (ir == NULL)
671  {
672  int order = Trans.OrderGrad(&el) + el.GetOrder();
673  ir = &IntRules.Get(el.GetGeomType(), order);
674  }
675 
676  Q.Eval(Q_nodal, Trans, el.GetNodes()); // sets the size of Q_nodal
677 
678  elmat = 0.0;
679  for (int i = 0; i < ir->GetNPoints(); i++)
680  {
681  const IntegrationPoint &ip = ir->IntPoint(i);
682  el.CalcDShape(ip, dshape);
683  el.CalcShape(ip, shape);
684 
685  Trans.SetIntPoint(&ip);
686  CalcAdjugate(Trans.Jacobian(), adjJ);
687 
688  Mult(dshape, adjJ, grad);
689 
690  double w = alpha * ip.weight;
691 
692  // elmat(k,l) += \sum_s w*shape(k)*Q_nodal(s,k)*grad(l,s)
693  for (int k = 0; k < nd; k++)
694  {
695  double wsk = w*shape(k);
696  for (int l = 0; l < nd; l++)
697  {
698  double a = 0.0;
699  for (int s = 0; s < dim; s++)
700  {
701  a += Q_nodal(s,k)*grad(l,s);
702  }
703  elmat(k,l) += wsk*a;
704  }
705  }
706  }
707 }
708 
709 
710 void VectorMassIntegrator::AssembleElementMatrix
712  DenseMatrix &elmat )
713 {
714  int nd = el.GetDof();
715  int spaceDim = Trans.GetSpaceDim();
716 
717  double norm;
718 
719  // Get vdim from VQ, MQ, or the space dimension
720  int vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : spaceDim);
721 
722  elmat.SetSize(nd*vdim);
723  shape.SetSize(nd);
724  partelmat.SetSize(nd);
725  if (VQ)
726  {
727  vec.SetSize(vdim);
728  }
729  else if (MQ)
730  {
731  mcoeff.SetSize(vdim);
732  }
733 
734  const IntegrationRule *ir = IntRule;
735  if (ir == NULL)
736  {
737  int order = 2 * el.GetOrder() + Trans.OrderW() + Q_order;
738 
739  if (el.Space() == FunctionSpace::rQk)
740  {
741  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
742  }
743  else
744  {
745  ir = &IntRules.Get(el.GetGeomType(), order);
746  }
747  }
748 
749  elmat = 0.0;
750  for (int s = 0; s < ir->GetNPoints(); s++)
751  {
752  const IntegrationPoint &ip = ir->IntPoint(s);
753  el.CalcShape(ip, shape);
754 
755  Trans.SetIntPoint (&ip);
756  norm = ip.weight * Trans.Weight();
757 
758  MultVVt(shape, partelmat);
759 
760  if (VQ)
761  {
762  VQ->Eval(vec, Trans, ip);
763  for (int k = 0; k < vdim; k++)
764  {
765  elmat.AddMatrix(norm*vec(k), partelmat, nd*k, nd*k);
766  }
767  }
768  else if (MQ)
769  {
770  MQ->Eval(mcoeff, Trans, ip);
771  for (int i = 0; i < vdim; i++)
772  for (int j = 0; j < vdim; j++)
773  {
774  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, nd*i, nd*j);
775  }
776  }
777  else
778  {
779  if (Q)
780  {
781  norm *= Q->Eval(Trans, ip);
782  }
783  partelmat *= norm;
784  for (int k = 0; k < vdim; k++)
785  {
786  elmat.AddMatrix(partelmat, nd*k, nd*k);
787  }
788  }
789  }
790 }
791 
792 void VectorMassIntegrator::AssembleElementMatrix2(
793  const FiniteElement &trial_fe, const FiniteElement &test_fe,
794  ElementTransformation &Trans, DenseMatrix &elmat)
795 {
796  int tr_nd = trial_fe.GetDof();
797  int te_nd = test_fe.GetDof();
798  int dim = trial_fe.GetDim();
799  int vdim;
800 
801  double norm;
802 
803  // Get vdim from the ElementTransformation Trans ?
804  vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : (dim));
805 
806  elmat.SetSize(te_nd*vdim, tr_nd*vdim);
807  shape.SetSize(tr_nd);
808  te_shape.SetSize(te_nd);
809  partelmat.SetSize(te_nd, tr_nd);
810  if (VQ)
811  {
812  vec.SetSize(vdim);
813  }
814  else if (MQ)
815  {
816  mcoeff.SetSize(vdim);
817  }
818 
819  const IntegrationRule *ir = IntRule;
820  if (ir == NULL)
821  {
822  int order = (trial_fe.GetOrder() + test_fe.GetOrder() +
823  Trans.OrderW() + Q_order);
824 
825  if (trial_fe.Space() == FunctionSpace::rQk)
826  {
827  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
828  }
829  else
830  {
831  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
832  }
833  }
834 
835  elmat = 0.0;
836  for (int s = 0; s < ir->GetNPoints(); s++)
837  {
838  const IntegrationPoint &ip = ir->IntPoint(s);
839  trial_fe.CalcShape(ip, shape);
840  test_fe.CalcShape(ip, te_shape);
841 
842  Trans.SetIntPoint(&ip);
843  norm = ip.weight * Trans.Weight();
844 
845  MultVWt(te_shape, shape, partelmat);
846 
847  if (VQ)
848  {
849  VQ->Eval(vec, Trans, ip);
850  for (int k = 0; k < vdim; k++)
851  {
852  elmat.AddMatrix(norm*vec(k), partelmat, te_nd*k, tr_nd*k);
853  }
854  }
855  else if (MQ)
856  {
857  MQ->Eval(mcoeff, Trans, ip);
858  for (int i = 0; i < vdim; i++)
859  for (int j = 0; j < vdim; j++)
860  {
861  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, te_nd*i, tr_nd*j);
862  }
863  }
864  else
865  {
866  if (Q)
867  {
868  norm *= Q->Eval(Trans, ip);
869  }
870  partelmat *= norm;
871  for (int k = 0; k < vdim; k++)
872  {
873  elmat.AddMatrix(partelmat, te_nd*k, tr_nd*k);
874  }
875  }
876  }
877 }
878 
879 void VectorFEDivergenceIntegrator::AssembleElementMatrix2(
880  const FiniteElement &trial_fe, const FiniteElement &test_fe,
881  ElementTransformation &Trans, DenseMatrix &elmat)
882 {
883  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
884 
885 #ifdef MFEM_THREAD_SAFE
886  Vector divshape(trial_nd), shape(test_nd);
887 #else
888  divshape.SetSize(trial_nd);
889  shape.SetSize(test_nd);
890 #endif
891 
892  elmat.SetSize(test_nd, trial_nd);
893 
894  const IntegrationRule *ir = IntRule;
895  if (ir == NULL)
896  {
897  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
898  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
899  }
900 
901  elmat = 0.0;
902  for (i = 0; i < ir->GetNPoints(); i++)
903  {
904  const IntegrationPoint &ip = ir->IntPoint(i);
905  trial_fe.CalcDivShape(ip, divshape);
906  test_fe.CalcShape(ip, shape);
907  double w = ip.weight;
908  if (Q)
909  {
910  Trans.SetIntPoint(&ip);
911  w *= Q->Eval(Trans, ip);
912  }
913  shape *= w;
914  AddMultVWt(shape, divshape, elmat);
915  }
916 }
917 
918 void VectorFECurlIntegrator::AssembleElementMatrix2(
919  const FiniteElement &trial_fe, const FiniteElement &test_fe,
920  ElementTransformation &Trans, DenseMatrix &elmat)
921 {
922  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
923  int dim = trial_fe.GetDim();
924 
925 #ifdef MFEM_THREAD_SAFE
926  DenseMatrix curlshapeTrial(trial_nd, dim);
927  DenseMatrix curlshapeTrial_dFT(trial_nd, dim);
928  DenseMatrix vshapeTest(test_nd, dim);
929 #else
930  curlshapeTrial.SetSize(trial_nd, dim);
931  curlshapeTrial_dFT.SetSize(trial_nd, dim);
932  vshapeTest.SetSize(test_nd, dim);
933 #endif
934 
935  elmat.SetSize(test_nd, trial_nd);
936 
937  const IntegrationRule *ir = IntRule;
938  if (ir == NULL)
939  {
940  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
941  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
942  }
943 
944  elmat = 0.0;
945  for (i = 0; i < ir->GetNPoints(); i++)
946  {
947  const IntegrationPoint &ip = ir->IntPoint(i);
948 
949  Trans.SetIntPoint(&ip);
950  trial_fe.CalcCurlShape(ip, curlshapeTrial);
951  MultABt(curlshapeTrial, Trans.Jacobian(), curlshapeTrial_dFT);
952  test_fe.CalcVShape(Trans, vshapeTest);
953  double w = ip.weight;
954 
955  if (Q)
956  {
957  w *= Q->Eval(Trans, ip);
958  }
959  vshapeTest *= w;
960  AddMultABt(vshapeTest, curlshapeTrial_dFT, elmat);
961  }
962 }
963 
964 void DerivativeIntegrator::AssembleElementMatrix2 (
965  const FiniteElement &trial_fe,
966  const FiniteElement &test_fe,
967  ElementTransformation &Trans,
968  DenseMatrix &elmat)
969 {
970  int dim = trial_fe.GetDim();
971  int trial_nd = trial_fe.GetDof();
972  int test_nd = test_fe.GetDof();
973 
974  int i, l;
975  double det;
976 
977  elmat.SetSize (test_nd,trial_nd);
978  dshape.SetSize (trial_nd,dim);
979  dshapedxt.SetSize(trial_nd,dim);
980  dshapedxi.SetSize(trial_nd);
981  invdfdx.SetSize(dim);
982  shape.SetSize (test_nd);
983 
984  const IntegrationRule *ir = IntRule;
985  if (ir == NULL)
986  {
987  int order;
988  if (trial_fe.Space() == FunctionSpace::Pk)
989  {
990  order = trial_fe.GetOrder() + test_fe.GetOrder() - 1;
991  }
992  else
993  {
994  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim;
995  }
996 
997  if (trial_fe.Space() == FunctionSpace::rQk)
998  {
999  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
1000  }
1001  else
1002  {
1003  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1004  }
1005  }
1006 
1007  elmat = 0.0;
1008  for (i = 0; i < ir->GetNPoints(); i++)
1009  {
1010  const IntegrationPoint &ip = ir->IntPoint(i);
1011 
1012  trial_fe.CalcDShape(ip, dshape);
1013 
1014  Trans.SetIntPoint (&ip);
1015  CalcInverse (Trans.Jacobian(), invdfdx);
1016  det = Trans.Weight();
1017  Mult (dshape, invdfdx, dshapedxt);
1018 
1019  test_fe.CalcShape(ip, shape);
1020 
1021  for (l = 0; l < trial_nd; l++)
1022  {
1023  dshapedxi(l) = dshapedxt(l,xi);
1024  }
1025 
1026  shape *= Q.Eval(Trans,ip) * det * ip.weight;
1027  AddMultVWt (shape, dshapedxi, elmat);
1028  }
1029 }
1030 
1031 void CurlCurlIntegrator::AssembleElementMatrix
1033  DenseMatrix &elmat )
1034 {
1035  int nd = el.GetDof();
1036  int dim = el.GetDim();
1037  int dimc = (dim == 3) ? 3 : 1;
1038  double w;
1039 
1040 #ifdef MFEM_THREAD_SAFE
1041  DenseMatrix curlshape(nd,dimc), curlshape_dFt(nd,dimc);
1042 #else
1043  curlshape.SetSize(nd,dimc);
1044  curlshape_dFt.SetSize(nd,dimc);
1045 #endif
1046  elmat.SetSize(nd);
1047 
1048  const IntegrationRule *ir = IntRule;
1049  if (ir == NULL)
1050  {
1051  int order;
1052  if (el.Space() == FunctionSpace::Pk)
1053  {
1054  order = 2*el.GetOrder() - 2;
1055  }
1056  else
1057  {
1058  order = 2*el.GetOrder();
1059  }
1060 
1061  ir = &IntRules.Get(el.GetGeomType(), order);
1062  }
1063 
1064  elmat = 0.0;
1065  for (int i = 0; i < ir->GetNPoints(); i++)
1066  {
1067  const IntegrationPoint &ip = ir->IntPoint(i);
1068 
1069  Trans.SetIntPoint (&ip);
1070 
1071  w = ip.weight / Trans.Weight();
1072 
1073  if ( dim == 3 )
1074  {
1075  el.CalcCurlShape(ip, curlshape);
1076  MultABt(curlshape, Trans.Jacobian(), curlshape_dFt);
1077  }
1078  else
1079  {
1080  el.CalcCurlShape(ip, curlshape_dFt);
1081  }
1082 
1083  if (Q)
1084  {
1085  w *= Q->Eval(Trans, ip);
1086  }
1087 
1088  AddMult_a_AAt(w, curlshape_dFt, elmat);
1089  }
1090 }
1091 
1092 void CurlCurlIntegrator
1093 ::ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans,
1094  Vector &u, const FiniteElement &fluxelem, Vector &flux,
1095  int with_coef)
1096 {
1097 #ifdef MFEM_THREAD_SAFE
1098  DenseMatrix projcurl;
1099 #endif
1100 
1101  fluxelem.ProjectCurl(el, Trans, projcurl);
1102 
1103  flux.SetSize(projcurl.Height());
1104  projcurl.Mult(u, flux);
1105 
1106  // TODO: Q, wcoef?
1107 }
1108 
1109 double CurlCurlIntegrator::ComputeFluxEnergy(const FiniteElement &fluxelem,
1110  ElementTransformation &Trans,
1111  Vector &flux, Vector *d_energy)
1112 {
1113  int nd = fluxelem.GetDof();
1114  int dim = fluxelem.GetDim();
1115 
1116 #ifdef MFEM_THREAD_SAFE
1117  DenseMatrix vshape;
1118 #endif
1119  vshape.SetSize(nd, dim);
1120  pointflux.SetSize(dim);
1121  if (d_energy) { vec.SetSize(dim); }
1122 
1123  int order = 2 * fluxelem.GetOrder(); // <--
1124  const IntegrationRule &ir = IntRules.Get(fluxelem.GetGeomType(), order);
1125 
1126  double energy = 0.0;
1127  if (d_energy) { *d_energy = 0.0; }
1128 
1129  Vector* pfluxes = NULL;
1130  if (d_energy)
1131  {
1132  pfluxes = new Vector[ir.GetNPoints()];
1133  }
1134 
1135  for (int i = 0; i < ir.GetNPoints(); i++)
1136  {
1137  const IntegrationPoint &ip = ir.IntPoint(i);
1138  Trans.SetIntPoint(&ip);
1139 
1140  fluxelem.CalcVShape(Trans, vshape);
1141  // fluxelem.CalcVShape(ip, vshape);
1142  vshape.MultTranspose(flux, pointflux);
1143 
1144  double w = Trans.Weight() * ip.weight;
1145 
1146  double e = w * (pointflux * pointflux);
1147 
1148  if (Q)
1149  {
1150  // TODO
1151  }
1152 
1153  energy += e;
1154 
1155 #if ANISO_EXPERIMENTAL
1156  if (d_energy)
1157  {
1158  pfluxes[i].SetSize(dim);
1159  Trans.Jacobian().MultTranspose(pointflux, pfluxes[i]);
1160 
1161  /*
1162  DenseMatrix Jadj(dim, dim);
1163  CalcAdjugate(Trans.Jacobian(), Jadj);
1164  pfluxes[i].SetSize(dim);
1165  Jadj.Mult(pointflux, pfluxes[i]);
1166  */
1167 
1168  // pfluxes[i] = pointflux;
1169  }
1170 #endif
1171  }
1172 
1173  if (d_energy)
1174  {
1175 #if ANISO_EXPERIMENTAL
1176  *d_energy = 0.0;
1177  Vector tmp;
1178 
1179  int n = (int) round(pow(ir.GetNPoints(), 1.0/3.0));
1180  MFEM_ASSERT(n*n*n == ir.GetNPoints(), "");
1181 
1182  // hack: get total variation of 'pointflux' in the x,y,z directions
1183  for (int k = 0; k < n; k++)
1184  for (int l = 0; l < n; l++)
1185  for (int m = 0; m < n; m++)
1186  {
1187  Vector &vec = pfluxes[(k*n + l)*n + m];
1188  if (m > 0)
1189  {
1190  tmp = vec; tmp -= pfluxes[(k*n + l)*n + (m-1)];
1191  (*d_energy)[0] += (tmp * tmp);
1192  }
1193  if (l > 0)
1194  {
1195  tmp = vec; tmp -= pfluxes[(k*n + (l-1))*n + m];
1196  (*d_energy)[1] += (tmp * tmp);
1197  }
1198  if (k > 0)
1199  {
1200  tmp = vec; tmp -= pfluxes[((k-1)*n + l)*n + m];
1201  (*d_energy)[2] += (tmp * tmp);
1202  }
1203  }
1204 #else
1205  *d_energy = 1.0;
1206 #endif
1207 
1208  delete [] pfluxes;
1209  }
1210 
1211  return energy;
1212 }
1213 
1214 void VectorCurlCurlIntegrator::AssembleElementMatrix(
1215  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1216 {
1217  int dim = el.GetDim();
1218  int dof = el.GetDof();
1219  int cld = (dim*(dim-1))/2;
1220 
1221 #ifdef MFEM_THREAD_SAFE
1222  DenseMatrix dshape_hat(dof, dim), dshape(dof, dim);
1223  DenseMatrix curlshape(dim*dof, cld), Jadj(dim);
1224 #else
1225  dshape_hat.SetSize(dof, dim);
1226  dshape.SetSize(dof, dim);
1227  curlshape.SetSize(dim*dof, cld);
1228  Jadj.SetSize(dim);
1229 #endif
1230 
1231  const IntegrationRule *ir = IntRule;
1232  if (ir == NULL)
1233  {
1234  // use the same integration rule as diffusion
1235  int order = 2 * Trans.OrderGrad(&el);
1236  ir = &IntRules.Get(el.GetGeomType(), order);
1237  }
1238 
1239  elmat = 0.0;
1240  for (int i = 0; i < ir->GetNPoints(); i++)
1241  {
1242  const IntegrationPoint &ip = ir->IntPoint(i);
1243  el.CalcDShape(ip, dshape_hat);
1244 
1245  Trans.SetIntPoint(&ip);
1246  CalcAdjugate(Trans.Jacobian(), Jadj);
1247  double w = ip.weight / Trans.Weight();
1248 
1249  Mult(dshape_hat, Jadj, dshape);
1250  dshape.GradToCurl(curlshape);
1251 
1252  if (Q)
1253  {
1254  w *= Q->Eval(Trans, ip);
1255  }
1256 
1257  AddMult_a_AAt(w, curlshape, elmat);
1258  }
1259 }
1260 
1261 double VectorCurlCurlIntegrator::GetElementEnergy(
1262  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
1263 {
1264  int dim = el.GetDim();
1265  int dof = el.GetDof();
1266 
1267 #ifdef MFEM_THREAD_SAFE
1268  DenseMatrix dshape_hat(dof, dim), Jadj(dim), grad_hat(dim), grad(dim);
1269 #else
1270  dshape_hat.SetSize(dof, dim);
1271 
1272  Jadj.SetSize(dim);
1273  grad_hat.SetSize(dim);
1274  grad.SetSize(dim);
1275 #endif
1276  DenseMatrix elfun_mat(elfun.GetData(), dof, dim);
1277 
1278  const IntegrationRule *ir = IntRule;
1279  if (ir == NULL)
1280  {
1281  // use the same integration rule as diffusion
1282  int order = 2 * Tr.OrderGrad(&el);
1283  ir = &IntRules.Get(el.GetGeomType(), order);
1284  }
1285 
1286  double energy = 0.;
1287  for (int i = 0; i < ir->GetNPoints(); i++)
1288  {
1289  const IntegrationPoint &ip = ir->IntPoint(i);
1290  el.CalcDShape(ip, dshape_hat);
1291 
1292  MultAtB(elfun_mat, dshape_hat, grad_hat);
1293 
1294  Tr.SetIntPoint(&ip);
1295  CalcAdjugate(Tr.Jacobian(), Jadj);
1296  double w = ip.weight / Tr.Weight();
1297 
1298  Mult(grad_hat, Jadj, grad);
1299 
1300  if (dim == 2)
1301  {
1302  double curl = grad(0,1) - grad(1,0);
1303  w *= curl * curl;
1304  }
1305  else
1306  {
1307  double curl_x = grad(2,1) - grad(1,2);
1308  double curl_y = grad(0,2) - grad(2,0);
1309  double curl_z = grad(1,0) - grad(0,1);
1310  w *= curl_x * curl_x + curl_y * curl_y + curl_z * curl_z;
1311  }
1312 
1313  if (Q)
1314  {
1315  w *= Q->Eval(Tr, ip);
1316  }
1317 
1318  energy += w;
1319  }
1320 
1321  elfun_mat.ClearExternalData();
1322 
1323  return 0.5 * energy;
1324 }
1325 
1326 
1327 void VectorFEMassIntegrator::AssembleElementMatrix(
1328  const FiniteElement &el,
1329  ElementTransformation &Trans,
1330  DenseMatrix &elmat)
1331 {
1332  int dof = el.GetDof();
1333  int spaceDim = Trans.GetSpaceDim();
1334 
1335  double w;
1336 
1337 #ifdef MFEM_THREAD_SAFE
1338  Vector D(VQ ? VQ->GetVDim() : 0);
1339  DenseMatrix trial_vshape(dof, spaceDim);
1340  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1341 #else
1342  trial_vshape.SetSize(dof, spaceDim);
1343  D.SetSize(VQ ? VQ->GetVDim() : 0);
1344  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1345 #endif
1346  DenseMatrix tmp(trial_vshape.Height(), K.Width());
1347 
1348  elmat.SetSize(dof);
1349  elmat = 0.0;
1350 
1351  const IntegrationRule *ir = IntRule;
1352  if (ir == NULL)
1353  {
1354  // int order = 2 * el.GetOrder();
1355  int order = Trans.OrderW() + 2 * el.GetOrder();
1356  ir = &IntRules.Get(el.GetGeomType(), order);
1357  }
1358 
1359  for (int i = 0; i < ir->GetNPoints(); i++)
1360  {
1361  const IntegrationPoint &ip = ir->IntPoint(i);
1362 
1363  Trans.SetIntPoint (&ip);
1364 
1365  el.CalcVShape(Trans, trial_vshape);
1366 
1367  w = ip.weight * Trans.Weight();
1368  if (MQ)
1369  {
1370  MQ->Eval(K, Trans, ip);
1371  K *= w;
1372  Mult(trial_vshape,K,tmp);
1373  AddMultABt(tmp,trial_vshape,elmat);
1374  }
1375  else if (VQ)
1376  {
1377  VQ->Eval(D, Trans, ip);
1378  D *= w;
1379  AddMultADAt(trial_vshape, D, elmat);
1380  }
1381  else
1382  {
1383  if (Q)
1384  {
1385  w *= Q -> Eval (Trans, ip);
1386  }
1387  AddMult_a_AAt (w, trial_vshape, elmat);
1388  }
1389  }
1390 }
1391 
1392 void VectorFEMassIntegrator::AssembleElementMatrix2(
1393  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1394  ElementTransformation &Trans, DenseMatrix &elmat)
1395 {
1396  if ( test_fe.GetRangeType() == FiniteElement::SCALAR && VQ )
1397  {
1398  // assume test_fe is scalar FE and trial_fe is vector FE
1399  int dim = test_fe.GetDim();
1400  int trial_dof = trial_fe.GetDof();
1401  int test_dof = test_fe.GetDof();
1402  double w;
1403 
1404  if (MQ)
1405  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1406  " is not implemented for tensor materials");
1407 
1408 #ifdef MFEM_THREAD_SAFE
1409  DenseMatrix trial_vshape(trial_dof, dim);
1410  Vector shape(test_dof);
1411  Vector D(dim);
1412 #else
1413  trial_vshape.SetSize(trial_dof, dim);
1414  shape.SetSize(test_dof);
1415  D.SetSize(dim);
1416 #endif
1417 
1418  elmat.SetSize (test_dof, trial_dof);
1419 
1420  const IntegrationRule *ir = IntRule;
1421  if (ir == NULL)
1422  {
1423  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1424  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1425  }
1426 
1427  elmat = 0.0;
1428  for (int i = 0; i < ir->GetNPoints(); i++)
1429  {
1430  const IntegrationPoint &ip = ir->IntPoint(i);
1431 
1432  Trans.SetIntPoint (&ip);
1433 
1434  trial_fe.CalcVShape(Trans, trial_vshape);
1435  test_fe.CalcShape(ip, shape);
1436 
1437  w = ip.weight * Trans.Weight();
1438  VQ->Eval(D, Trans, ip);
1439  D *= w;
1440 
1441  for (int d = 0; d < dim; d++)
1442  {
1443  for (int j = 0; j < test_dof; j++)
1444  {
1445  for (int k = 0; k < trial_dof; k++)
1446  {
1447  elmat(j, k) += D[d] * shape(j) * trial_vshape(k, d);
1448  }
1449  }
1450  }
1451  }
1452  }
1453  else if ( test_fe.GetRangeType() == FiniteElement::SCALAR )
1454  {
1455  // assume test_fe is scalar FE and trial_fe is vector FE
1456  int dim = test_fe.GetDim();
1457  int trial_dof = trial_fe.GetDof();
1458  int test_dof = test_fe.GetDof();
1459  double w;
1460 
1461  if (VQ || MQ)
1462  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1463  " is not implemented for vector/tensor permeability");
1464 
1465 #ifdef MFEM_THREAD_SAFE
1466  DenseMatrix trial_vshape(trial_dof, dim);
1467  Vector shape(test_dof);
1468 #else
1469  trial_vshape.SetSize(trial_dof, dim);
1470  shape.SetSize(test_dof);
1471 #endif
1472 
1473  elmat.SetSize (dim*test_dof, trial_dof);
1474 
1475  const IntegrationRule *ir = IntRule;
1476  if (ir == NULL)
1477  {
1478  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1479  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1480  }
1481 
1482  elmat = 0.0;
1483  for (int i = 0; i < ir->GetNPoints(); i++)
1484  {
1485  const IntegrationPoint &ip = ir->IntPoint(i);
1486 
1487  Trans.SetIntPoint (&ip);
1488 
1489  trial_fe.CalcVShape(Trans, trial_vshape);
1490  test_fe.CalcShape(ip, shape);
1491 
1492  w = ip.weight * Trans.Weight();
1493  if (Q)
1494  {
1495  w *= Q -> Eval (Trans, ip);
1496  }
1497 
1498  for (int d = 0; d < dim; d++)
1499  {
1500  for (int j = 0; j < test_dof; j++)
1501  {
1502  for (int k = 0; k < trial_dof; k++)
1503  {
1504  elmat(d * test_dof + j, k) += w * shape(j) * trial_vshape(k, d);
1505  }
1506  }
1507  }
1508  }
1509  }
1510  else
1511  {
1512  // assume both test_fe and trial_fe are vector FE
1513  int dim = test_fe.GetDim();
1514  int trial_dof = trial_fe.GetDof();
1515  int test_dof = test_fe.GetDof();
1516  double w;
1517 
1518  if (VQ || MQ)
1519  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1520  " is not implemented for vector/tensor permeability");
1521 
1522 #ifdef MFEM_THREAD_SAFE
1523  DenseMatrix trial_vshape(trial_dof, dim);
1524  DenseMatrix test_vshape(test_dof,dim);
1525 #else
1526  trial_vshape.SetSize(trial_dof, dim);
1527  test_vshape.SetSize(test_dof,dim);
1528 #endif
1529 
1530  elmat.SetSize (test_dof, trial_dof);
1531 
1532  const IntegrationRule *ir = IntRule;
1533  if (ir == NULL)
1534  {
1535  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1536  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1537  }
1538 
1539  elmat = 0.0;
1540  for (int i = 0; i < ir->GetNPoints(); i++)
1541  {
1542  const IntegrationPoint &ip = ir->IntPoint(i);
1543 
1544  Trans.SetIntPoint (&ip);
1545 
1546  trial_fe.CalcVShape(Trans, trial_vshape);
1547  test_fe.CalcVShape(Trans, test_vshape);
1548 
1549  w = ip.weight * Trans.Weight();
1550  if (Q)
1551  {
1552  w *= Q -> Eval (Trans, ip);
1553  }
1554 
1555  for (int d = 0; d < dim; d++)
1556  {
1557  for (int j = 0; j < test_dof; j++)
1558  {
1559  for (int k = 0; k < trial_dof; k++)
1560  {
1561  elmat(j, k) += w * test_vshape(j, d) * trial_vshape(k, d);
1562  }
1563  }
1564  }
1565  }
1566  }
1567 }
1568 
1569 void VectorDivergenceIntegrator::AssembleElementMatrix2(
1570  const FiniteElement &trial_fe,
1571  const FiniteElement &test_fe,
1572  ElementTransformation &Trans,
1573  DenseMatrix &elmat)
1574 {
1575  int dim = trial_fe.GetDim();
1576  int trial_dof = trial_fe.GetDof();
1577  int test_dof = test_fe.GetDof();
1578  double c;
1579 
1580  dshape.SetSize (trial_dof, dim);
1581  gshape.SetSize (trial_dof, dim);
1582  Jadj.SetSize (dim);
1583  divshape.SetSize (dim*trial_dof);
1584  shape.SetSize (test_dof);
1585 
1586  elmat.SetSize (test_dof, dim*trial_dof);
1587 
1588  const IntegrationRule *ir = IntRule;
1589  if (ir == NULL)
1590  {
1591  int order = Trans.OrderGrad(&trial_fe) + test_fe.GetOrder();
1592  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1593  }
1594 
1595  elmat = 0.0;
1596 
1597  for (int i = 0; i < ir -> GetNPoints(); i++)
1598  {
1599  const IntegrationPoint &ip = ir->IntPoint(i);
1600 
1601  trial_fe.CalcDShape (ip, dshape);
1602  test_fe.CalcShape (ip, shape);
1603 
1604  Trans.SetIntPoint (&ip);
1605  CalcAdjugate(Trans.Jacobian(), Jadj);
1606 
1607  Mult (dshape, Jadj, gshape);
1608 
1609  gshape.GradToDiv (divshape);
1610 
1611  c = ip.weight;
1612  if (Q)
1613  {
1614  c *= Q -> Eval (Trans, ip);
1615  }
1616 
1617  // elmat += c * shape * divshape ^ t
1618  shape *= c;
1619  AddMultVWt (shape, divshape, elmat);
1620  }
1621 }
1622 
1623 
1624 void DivDivIntegrator::AssembleElementMatrix(
1625  const FiniteElement &el,
1626  ElementTransformation &Trans,
1627  DenseMatrix &elmat)
1628 {
1629  int dof = el.GetDof();
1630  double c;
1631 
1632 #ifdef MFEM_THREAD_SAFE
1633  Vector divshape(dof);
1634 #else
1635  divshape.SetSize(dof);
1636 #endif
1637  elmat.SetSize(dof);
1638 
1639  const IntegrationRule *ir = IntRule;
1640  if (ir == NULL)
1641  {
1642  int order = 2 * el.GetOrder() - 2; // <--- OK for RTk
1643  ir = &IntRules.Get(el.GetGeomType(), order);
1644  }
1645 
1646  elmat = 0.0;
1647 
1648  for (int i = 0; i < ir -> GetNPoints(); i++)
1649  {
1650  const IntegrationPoint &ip = ir->IntPoint(i);
1651 
1652  el.CalcDivShape (ip, divshape);
1653 
1654  Trans.SetIntPoint (&ip);
1655  c = ip.weight / Trans.Weight();
1656 
1657  if (Q)
1658  {
1659  c *= Q -> Eval (Trans, ip);
1660  }
1661 
1662  // elmat += c * divshape * divshape ^ t
1663  AddMult_a_VVt (c, divshape, elmat);
1664  }
1665 }
1666 
1667 
1668 void VectorDiffusionIntegrator::AssembleElementMatrix(
1669  const FiniteElement &el,
1670  ElementTransformation &Trans,
1671  DenseMatrix &elmat)
1672 {
1673  int dim = el.GetDim();
1674  int dof = el.GetDof();
1675 
1676  double norm;
1677 
1678  elmat.SetSize (dim * dof);
1679 
1680  Jinv. SetSize (dim);
1681  dshape.SetSize (dof, dim);
1682  gshape.SetSize (dof, dim);
1683  pelmat.SetSize (dof);
1684 
1685  const IntegrationRule *ir = IntRule;
1686  if (ir == NULL)
1687  {
1688  // integrand is rational function if det(J) is not constant
1689  int order = 2 * Trans.OrderGrad(&el); // order of the numerator
1690  if (el.Space() == FunctionSpace::rQk)
1691  {
1692  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
1693  }
1694  else
1695  {
1696  ir = &IntRules.Get(el.GetGeomType(), order);
1697  }
1698  }
1699 
1700  elmat = 0.0;
1701 
1702  for (int i = 0; i < ir -> GetNPoints(); i++)
1703  {
1704  const IntegrationPoint &ip = ir->IntPoint(i);
1705 
1706  el.CalcDShape (ip, dshape);
1707 
1708  Trans.SetIntPoint (&ip);
1709  norm = ip.weight * Trans.Weight();
1710  CalcInverse (Trans.Jacobian(), Jinv);
1711 
1712  Mult (dshape, Jinv, gshape);
1713 
1714  MultAAt (gshape, pelmat);
1715 
1716  if (Q)
1717  {
1718  norm *= Q -> Eval (Trans, ip);
1719  }
1720 
1721  pelmat *= norm;
1722 
1723  for (int d = 0; d < dim; d++)
1724  {
1725  for (int k = 0; k < dof; k++)
1726  for (int l = 0; l < dof; l++)
1727  {
1728  elmat (dof*d+k, dof*d+l) += pelmat (k, l);
1729  }
1730  }
1731  }
1732 }
1733 
1734 
1735 void ElasticityIntegrator::AssembleElementMatrix(
1736  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1737 {
1738  int dof = el.GetDof();
1739  int dim = el.GetDim();
1740  double w, L, M;
1741 
1742 #ifdef MFEM_THREAD_SAFE
1743  DenseMatrix dshape(dof, dim), Jinv(dim), gshape(dof, dim), pelmat(dof);
1744  Vector divshape(dim*dof);
1745 #else
1746  Jinv.SetSize(dim);
1747  dshape.SetSize(dof, dim);
1748  gshape.SetSize(dof, dim);
1749  pelmat.SetSize(dof);
1750  divshape.SetSize(dim*dof);
1751 #endif
1752 
1753  elmat.SetSize(dof * dim);
1754 
1755  const IntegrationRule *ir = IntRule;
1756  if (ir == NULL)
1757  {
1758  int order = 2 * Trans.OrderGrad(&el); // correct order?
1759  ir = &IntRules.Get(el.GetGeomType(), order);
1760  }
1761 
1762  elmat = 0.0;
1763 
1764  for (int i = 0; i < ir -> GetNPoints(); i++)
1765  {
1766  const IntegrationPoint &ip = ir->IntPoint(i);
1767 
1768  el.CalcDShape(ip, dshape);
1769 
1770  Trans.SetIntPoint(&ip);
1771  w = ip.weight * Trans.Weight();
1772  CalcInverse(Trans.Jacobian(), Jinv);
1773  Mult(dshape, Jinv, gshape);
1774  MultAAt(gshape, pelmat);
1775  gshape.GradToDiv (divshape);
1776 
1777  M = mu->Eval(Trans, ip);
1778  if (lambda)
1779  {
1780  L = lambda->Eval(Trans, ip);
1781  }
1782  else
1783  {
1784  L = q_lambda * M;
1785  M = q_mu * M;
1786  }
1787 
1788  if (L != 0.0)
1789  {
1790  AddMult_a_VVt(L * w, divshape, elmat);
1791  }
1792 
1793  if (M != 0.0)
1794  {
1795  for (int d = 0; d < dim; d++)
1796  {
1797  for (int k = 0; k < dof; k++)
1798  for (int l = 0; l < dof; l++)
1799  {
1800  elmat (dof*d+k, dof*d+l) += (M * w) * pelmat(k, l);
1801  }
1802  }
1803  for (int i = 0; i < dim; i++)
1804  for (int j = 0; j < dim; j++)
1805  {
1806  for (int k = 0; k < dof; k++)
1807  for (int l = 0; l < dof; l++)
1808  elmat(dof*i+k, dof*j+l) +=
1809  (M * w) * gshape(k, j) * gshape(l, i);
1810  // + (L * w) * gshape(k, i) * gshape(l, j)
1811  }
1812  }
1813  }
1814 }
1815 
1816 void DGTraceIntegrator::AssembleFaceMatrix(const FiniteElement &el1,
1817  const FiniteElement &el2,
1819  DenseMatrix &elmat)
1820 {
1821  int dim, ndof1, ndof2;
1822 
1823  double un, a, b, w;
1824 
1825  dim = el1.GetDim();
1826  ndof1 = el1.GetDof();
1827  Vector vu(dim), nor(dim);
1828 
1829  if (Trans.Elem2No >= 0)
1830  {
1831  ndof2 = el2.GetDof();
1832  }
1833  else
1834  {
1835  ndof2 = 0;
1836  }
1837 
1838  shape1.SetSize(ndof1);
1839  shape2.SetSize(ndof2);
1840  elmat.SetSize(ndof1 + ndof2);
1841  elmat = 0.0;
1842 
1843  const IntegrationRule *ir = IntRule;
1844  if (ir == NULL)
1845  {
1846  int order;
1847  // Assuming order(u)==order(mesh)
1848  if (Trans.Elem2No >= 0)
1849  order = (min(Trans.Elem1->OrderW(), Trans.Elem2->OrderW()) +
1850  2*max(el1.GetOrder(), el2.GetOrder()));
1851  else
1852  {
1853  order = Trans.Elem1->OrderW() + 2*el1.GetOrder();
1854  }
1855  if (el1.Space() == FunctionSpace::Pk)
1856  {
1857  order++;
1858  }
1859  ir = &IntRules.Get(Trans.FaceGeom, order);
1860  }
1861 
1862  for (int p = 0; p < ir->GetNPoints(); p++)
1863  {
1864  const IntegrationPoint &ip = ir->IntPoint(p);
1865  IntegrationPoint eip1, eip2;
1866  Trans.Loc1.Transform(ip, eip1);
1867  if (ndof2)
1868  {
1869  Trans.Loc2.Transform(ip, eip2);
1870  }
1871  el1.CalcShape(eip1, shape1);
1872 
1873  Trans.Face->SetIntPoint(&ip);
1874  Trans.Elem1->SetIntPoint(&eip1);
1875 
1876  u->Eval(vu, *Trans.Elem1, eip1);
1877 
1878  if (dim == 1)
1879  {
1880  nor(0) = 2*eip1.x - 1.0;
1881  }
1882  else
1883  {
1884  CalcOrtho(Trans.Face->Jacobian(), nor);
1885  }
1886 
1887  un = vu * nor;
1888  a = 0.5 * alpha * un;
1889  b = beta * fabs(un);
1890  // note: if |alpha/2|==|beta| then |a|==|b|, i.e. (a==b) or (a==-b)
1891  // and therefore two blocks in the element matrix contribution
1892  // (from the current quadrature point) are 0
1893 
1894  if (rho)
1895  {
1896  double rho_p;
1897  if (un >= 0.0 && ndof2)
1898  {
1899  Trans.Elem2->SetIntPoint(&eip2);
1900  rho_p = rho->Eval(*Trans.Elem2, eip2);
1901  }
1902  else
1903  {
1904  rho_p = rho->Eval(*Trans.Elem1, eip1);
1905  }
1906  a *= rho_p;
1907  b *= rho_p;
1908  }
1909 
1910  w = ip.weight * (a+b);
1911  if (w != 0.0)
1912  {
1913  for (int i = 0; i < ndof1; i++)
1914  for (int j = 0; j < ndof1; j++)
1915  {
1916  elmat(i, j) += w * shape1(i) * shape1(j);
1917  }
1918  }
1919 
1920  if (ndof2)
1921  {
1922  el2.CalcShape(eip2, shape2);
1923 
1924  if (w != 0.0)
1925  for (int i = 0; i < ndof2; i++)
1926  for (int j = 0; j < ndof1; j++)
1927  {
1928  elmat(ndof1+i, j) -= w * shape2(i) * shape1(j);
1929  }
1930 
1931  w = ip.weight * (b-a);
1932  if (w != 0.0)
1933  {
1934  for (int i = 0; i < ndof2; i++)
1935  for (int j = 0; j < ndof2; j++)
1936  {
1937  elmat(ndof1+i, ndof1+j) += w * shape2(i) * shape2(j);
1938  }
1939 
1940  for (int i = 0; i < ndof1; i++)
1941  for (int j = 0; j < ndof2; j++)
1942  {
1943  elmat(i, ndof1+j) -= w * shape1(i) * shape2(j);
1944  }
1945  }
1946  }
1947  }
1948 }
1949 
1950 void DGDiffusionIntegrator::AssembleFaceMatrix(
1951  const FiniteElement &el1, const FiniteElement &el2,
1952  FaceElementTransformations &Trans, DenseMatrix &elmat)
1953 {
1954  int dim, ndof1, ndof2, ndofs;
1955  bool kappa_is_nonzero = (kappa != 0.);
1956  double w, wq = 0.0;
1957 
1958  dim = el1.GetDim();
1959  ndof1 = el1.GetDof();
1960 
1961  nor.SetSize(dim);
1962  nh.SetSize(dim);
1963  ni.SetSize(dim);
1964  adjJ.SetSize(dim);
1965  if (MQ)
1966  {
1967  mq.SetSize(dim);
1968  }
1969 
1970  shape1.SetSize(ndof1);
1971  dshape1.SetSize(ndof1, dim);
1972  dshape1dn.SetSize(ndof1);
1973  if (Trans.Elem2No >= 0)
1974  {
1975  ndof2 = el2.GetDof();
1976  shape2.SetSize(ndof2);
1977  dshape2.SetSize(ndof2, dim);
1978  dshape2dn.SetSize(ndof2);
1979  }
1980  else
1981  {
1982  ndof2 = 0;
1983  }
1984 
1985  ndofs = ndof1 + ndof2;
1986  elmat.SetSize(ndofs);
1987  elmat = 0.0;
1988  if (kappa_is_nonzero)
1989  {
1990  jmat.SetSize(ndofs);
1991  jmat = 0.;
1992  }
1993 
1994  const IntegrationRule *ir = IntRule;
1995  if (ir == NULL)
1996  {
1997  // a simple choice for the integration order; is this OK?
1998  int order;
1999  if (ndof2)
2000  {
2001  order = 2*max(el1.GetOrder(), el2.GetOrder());
2002  }
2003  else
2004  {
2005  order = 2*el1.GetOrder();
2006  }
2007  ir = &IntRules.Get(Trans.FaceGeom, order);
2008  }
2009 
2010  // assemble: < {(Q \nabla u).n},[v] > --> elmat
2011  // kappa < {h^{-1} Q} [u],[v] > --> jmat
2012  for (int p = 0; p < ir->GetNPoints(); p++)
2013  {
2014  const IntegrationPoint &ip = ir->IntPoint(p);
2015  IntegrationPoint eip1, eip2;
2016 
2017  Trans.Loc1.Transform(ip, eip1);
2018  Trans.Face->SetIntPoint(&ip);
2019  if (dim == 1)
2020  {
2021  nor(0) = 2*eip1.x - 1.0;
2022  }
2023  else
2024  {
2025  CalcOrtho(Trans.Face->Jacobian(), nor);
2026  }
2027 
2028  el1.CalcShape(eip1, shape1);
2029  el1.CalcDShape(eip1, dshape1);
2030  Trans.Elem1->SetIntPoint(&eip1);
2031  w = ip.weight/Trans.Elem1->Weight();
2032  if (ndof2)
2033  {
2034  w /= 2;
2035  }
2036  if (!MQ)
2037  {
2038  if (Q)
2039  {
2040  w *= Q->Eval(*Trans.Elem1, eip1);
2041  }
2042  ni.Set(w, nor);
2043  }
2044  else
2045  {
2046  nh.Set(w, nor);
2047  MQ->Eval(mq, *Trans.Elem1, eip1);
2048  mq.MultTranspose(nh, ni);
2049  }
2050  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
2051  adjJ.Mult(ni, nh);
2052  if (kappa_is_nonzero)
2053  {
2054  wq = ni * nor;
2055  }
2056  // Note: in the jump term, we use 1/h1 = |nor|/det(J1) which is
2057  // independent of Loc1 and always gives the size of element 1 in
2058  // direction perpendicular to the face. Indeed, for linear transformation
2059  // |nor|=measure(face)/measure(ref. face),
2060  // det(J1)=measure(element)/measure(ref. element),
2061  // and the ratios measure(ref. element)/measure(ref. face) are
2062  // compatible for all element/face pairs.
2063  // For example: meas(ref. tetrahedron)/meas(ref. triangle) = 1/3, and
2064  // for any tetrahedron vol(tet)=(1/3)*height*area(base).
2065  // For interior faces: q_e/h_e=(q1/h1+q2/h2)/2.
2066 
2067  dshape1.Mult(nh, dshape1dn);
2068  for (int i = 0; i < ndof1; i++)
2069  for (int j = 0; j < ndof1; j++)
2070  {
2071  elmat(i, j) += shape1(i) * dshape1dn(j);
2072  }
2073 
2074  if (ndof2)
2075  {
2076  Trans.Loc2.Transform(ip, eip2);
2077  el2.CalcShape(eip2, shape2);
2078  el2.CalcDShape(eip2, dshape2);
2079  Trans.Elem2->SetIntPoint(&eip2);
2080  w = ip.weight/2/Trans.Elem2->Weight();
2081  if (!MQ)
2082  {
2083  if (Q)
2084  {
2085  w *= Q->Eval(*Trans.Elem2, eip2);
2086  }
2087  ni.Set(w, nor);
2088  }
2089  else
2090  {
2091  nh.Set(w, nor);
2092  MQ->Eval(mq, *Trans.Elem2, eip2);
2093  mq.MultTranspose(nh, ni);
2094  }
2095  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
2096  adjJ.Mult(ni, nh);
2097  if (kappa_is_nonzero)
2098  {
2099  wq += ni * nor;
2100  }
2101 
2102  dshape2.Mult(nh, dshape2dn);
2103 
2104  for (int i = 0; i < ndof1; i++)
2105  for (int j = 0; j < ndof2; j++)
2106  {
2107  elmat(i, ndof1 + j) += shape1(i) * dshape2dn(j);
2108  }
2109 
2110  for (int i = 0; i < ndof2; i++)
2111  for (int j = 0; j < ndof1; j++)
2112  {
2113  elmat(ndof1 + i, j) -= shape2(i) * dshape1dn(j);
2114  }
2115 
2116  for (int i = 0; i < ndof2; i++)
2117  for (int j = 0; j < ndof2; j++)
2118  {
2119  elmat(ndof1 + i, ndof1 + j) -= shape2(i) * dshape2dn(j);
2120  }
2121  }
2122 
2123  if (kappa_is_nonzero)
2124  {
2125  // only assemble the lower triangular part of jmat
2126  wq *= kappa;
2127  for (int i = 0; i < ndof1; i++)
2128  {
2129  const double wsi = wq*shape1(i);
2130  for (int j = 0; j <= i; j++)
2131  {
2132  jmat(i, j) += wsi * shape1(j);
2133  }
2134  }
2135  if (ndof2)
2136  {
2137  for (int i = 0; i < ndof2; i++)
2138  {
2139  const int i2 = ndof1 + i;
2140  const double wsi = wq*shape2(i);
2141  for (int j = 0; j < ndof1; j++)
2142  {
2143  jmat(i2, j) -= wsi * shape1(j);
2144  }
2145  for (int j = 0; j <= i; j++)
2146  {
2147  jmat(i2, ndof1 + j) += wsi * shape2(j);
2148  }
2149  }
2150  }
2151  }
2152  }
2153 
2154  // elmat := -elmat + sigma*elmat^t + jmat
2155  if (kappa_is_nonzero)
2156  {
2157  for (int i = 0; i < ndofs; i++)
2158  {
2159  for (int j = 0; j < i; j++)
2160  {
2161  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
2162  elmat(i,j) = sigma*aji - aij + mij;
2163  elmat(j,i) = sigma*aij - aji + mij;
2164  }
2165  elmat(i,i) = (sigma - 1.)*elmat(i,i) + jmat(i,i);
2166  }
2167  }
2168  else
2169  {
2170  for (int i = 0; i < ndofs; i++)
2171  {
2172  for (int j = 0; j < i; j++)
2173  {
2174  double aij = elmat(i,j), aji = elmat(j,i);
2175  elmat(i,j) = sigma*aji - aij;
2176  elmat(j,i) = sigma*aij - aji;
2177  }
2178  elmat(i,i) *= (sigma - 1.);
2179  }
2180  }
2181 }
2182 
2183 void TraceJumpIntegrator::AssembleFaceMatrix(
2184  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2185  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2186  DenseMatrix &elmat)
2187 {
2188  int i, j, face_ndof, ndof1, ndof2;
2189  int order;
2190 
2191  double w;
2192 
2193  face_ndof = trial_face_fe.GetDof();
2194  ndof1 = test_fe1.GetDof();
2195 
2196  face_shape.SetSize(face_ndof);
2197  shape1.SetSize(ndof1);
2198 
2199  if (Trans.Elem2No >= 0)
2200  {
2201  ndof2 = test_fe2.GetDof();
2202  shape2.SetSize(ndof2);
2203  }
2204  else
2205  {
2206  ndof2 = 0;
2207  }
2208 
2209  elmat.SetSize(ndof1 + ndof2, face_ndof);
2210  elmat = 0.0;
2211 
2212  const IntegrationRule *ir = IntRule;
2213  if (ir == NULL)
2214  {
2215  if (Trans.Elem2No >= 0)
2216  {
2217  order = max(test_fe1.GetOrder(), test_fe2.GetOrder());
2218  }
2219  else
2220  {
2221  order = test_fe1.GetOrder();
2222  }
2223  order += trial_face_fe.GetOrder();
2224  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2225  {
2226  order += Trans.Face->OrderW();
2227  }
2228  ir = &IntRules.Get(Trans.FaceGeom, order);
2229  }
2230 
2231  for (int p = 0; p < ir->GetNPoints(); p++)
2232  {
2233  const IntegrationPoint &ip = ir->IntPoint(p);
2234  IntegrationPoint eip1, eip2;
2235  // Trace finite element shape function
2236  Trans.Face->SetIntPoint(&ip);
2237  trial_face_fe.CalcShape(ip, face_shape);
2238  // Side 1 finite element shape function
2239  Trans.Loc1.Transform(ip, eip1);
2240  test_fe1.CalcShape(eip1, shape1);
2241  Trans.Elem1->SetIntPoint(&eip1);
2242  if (ndof2)
2243  {
2244  // Side 2 finite element shape function
2245  Trans.Loc2.Transform(ip, eip2);
2246  test_fe2.CalcShape(eip2, shape2);
2247  Trans.Elem2->SetIntPoint(&eip2);
2248  }
2249  w = ip.weight;
2250  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2251  {
2252  w *= Trans.Face->Weight();
2253  }
2254  face_shape *= w;
2255  for (i = 0; i < ndof1; i++)
2256  for (j = 0; j < face_ndof; j++)
2257  {
2258  elmat(i, j) += shape1(i) * face_shape(j);
2259  }
2260  if (ndof2)
2261  {
2262  // Subtract contribution from side 2
2263  for (i = 0; i < ndof2; i++)
2264  for (j = 0; j < face_ndof; j++)
2265  {
2266  elmat(ndof1+i, j) -= shape2(i) * face_shape(j);
2267  }
2268  }
2269  }
2270 }
2271 
2272 void NormalTraceJumpIntegrator::AssembleFaceMatrix(
2273  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2274  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2275  DenseMatrix &elmat)
2276 {
2277  int i, j, face_ndof, ndof1, ndof2, dim;
2278  int order;
2279 
2280  MFEM_VERIFY(trial_face_fe.GetMapType() == FiniteElement::VALUE, "");
2281 
2282  face_ndof = trial_face_fe.GetDof();
2283  ndof1 = test_fe1.GetDof();
2284  dim = test_fe1.GetDim();
2285 
2286  face_shape.SetSize(face_ndof);
2287  normal.SetSize(dim);
2288  shape1.SetSize(ndof1,dim);
2289  shape1_n.SetSize(ndof1);
2290 
2291  if (Trans.Elem2No >= 0)
2292  {
2293  ndof2 = test_fe2.GetDof();
2294  shape2.SetSize(ndof2,dim);
2295  shape2_n.SetSize(ndof2);
2296  }
2297  else
2298  {
2299  ndof2 = 0;
2300  }
2301 
2302  elmat.SetSize(ndof1 + ndof2, face_ndof);
2303  elmat = 0.0;
2304 
2305  const IntegrationRule *ir = IntRule;
2306  if (ir == NULL)
2307  {
2308  if (Trans.Elem2No >= 0)
2309  {
2310  order = max(test_fe1.GetOrder(), test_fe2.GetOrder()) - 1;
2311  }
2312  else
2313  {
2314  order = test_fe1.GetOrder() - 1;
2315  }
2316  order += trial_face_fe.GetOrder();
2317  ir = &IntRules.Get(Trans.FaceGeom, order);
2318  }
2319 
2320  for (int p = 0; p < ir->GetNPoints(); p++)
2321  {
2322  const IntegrationPoint &ip = ir->IntPoint(p);
2323  IntegrationPoint eip1, eip2;
2324  // Trace finite element shape function
2325  trial_face_fe.CalcShape(ip, face_shape);
2326  Trans.Loc1.Transf.SetIntPoint(&ip);
2327  CalcOrtho(Trans.Loc1.Transf.Jacobian(), normal);
2328  // Side 1 finite element shape function
2329  Trans.Loc1.Transform(ip, eip1);
2330  test_fe1.CalcVShape(eip1, shape1);
2331  shape1.Mult(normal, shape1_n);
2332  if (ndof2)
2333  {
2334  // Side 2 finite element shape function
2335  Trans.Loc2.Transform(ip, eip2);
2336  test_fe2.CalcVShape(eip2, shape2);
2337  Trans.Loc2.Transf.SetIntPoint(&ip);
2338  CalcOrtho(Trans.Loc2.Transf.Jacobian(), normal);
2339  shape2.Mult(normal, shape2_n);
2340  }
2341  face_shape *= ip.weight;
2342  for (i = 0; i < ndof1; i++)
2343  for (j = 0; j < face_ndof; j++)
2344  {
2345  elmat(i, j) -= shape1_n(i) * face_shape(j);
2346  }
2347  if (ndof2)
2348  {
2349  // Subtract contribution from side 2
2350  for (i = 0; i < ndof2; i++)
2351  for (j = 0; j < face_ndof; j++)
2352  {
2353  elmat(ndof1+i, j) += shape2_n(i) * face_shape(j);
2354  }
2355  }
2356  }
2357 }
2358 
2359 
2360 void NormalInterpolator::AssembleElementMatrix2(
2361  const FiniteElement &dom_fe, const FiniteElement &ran_fe,
2362  ElementTransformation &Trans, DenseMatrix &elmat)
2363 {
2364  int spaceDim = Trans.GetSpaceDim();
2365  elmat.SetSize(ran_fe.GetDof(), spaceDim*dom_fe.GetDof());
2366  Vector n(spaceDim), shape(dom_fe.GetDof());
2367 
2368  const IntegrationRule &ran_nodes = ran_fe.GetNodes();
2369  for (int i = 0; i < ran_nodes.Size(); i++)
2370  {
2371  const IntegrationPoint &ip = ran_nodes.IntPoint(i);
2372  Trans.SetIntPoint(&ip);
2373  CalcOrtho(Trans.Jacobian(), n);
2374  dom_fe.CalcShape(ip, shape);
2375  for (int j = 0; j < shape.Size(); j++)
2376  {
2377  for (int d = 0; d < spaceDim; d++)
2378  {
2379  elmat(i, j+d*shape.Size()) = shape(j)*n(d);
2380  }
2381  }
2382  }
2383 }
2384 
2385 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:228
Abstract class for Finite Elements.
Definition: fe.hpp:44
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:3218
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:82
ElementTransformation * Face
Definition: eltrans.hpp:123
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3525
Class for integration rule.
Definition: intrules.hpp:83
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:30
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:3503
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:235
void SetSize(int s)
Resize the vector if the new size is different.
Definition: vector.hpp:263
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:451
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:271
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:35
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2904
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:91
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:124
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:94
int GetMapType() const
Definition: fe.hpp:98
ElementTransformation * Elem2
Definition: eltrans.hpp:123
double * GetData() const
Definition: vector.hpp:90
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3134
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:108
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:231
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:193
int dim
Definition: ex3.cpp:47
IntegrationRules IntRules(0)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:292
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:124
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:566
const IntegrationRule & GetNodes() const
Definition: fe.hpp:113
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:3568
virtual double Weight()=0
virtual const DenseMatrix & Jacobian()
Definition: eltrans.cpp:52
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:3012
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:85
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:2294
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:3494
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:3341
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
void mfem_error(const char *msg)
Definition: error.cpp:23
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:243
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:3160
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:44
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:2559
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:51
Class for integration point with weight.
Definition: intrules.hpp:25
virtual int OrderGrad(const FiniteElement *fe)=0
order of adj(J)^t.grad(fi)
IsoparametricTransformation Transf
Definition: eltrans.hpp:114
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:2358
ElementTransformation * Elem1
Definition: eltrans.hpp:123
const double alpha
Definition: ex15.cpp:337
double kappa
Definition: ex3.cpp:46
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:3398
Vector data type.
Definition: vector.hpp:33
IntegrationRules RefinedIntRules(1)
A global object with all refined integration rules.
Definition: intrules.hpp:295
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
Definition: densemat.cpp:3174
virtual int GetSpaceDim()=0
virtual const DenseMatrix & Jacobian()=0
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:71
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
int GetRangeType() const
Definition: fe.hpp:96
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:3455