MFEM  v3.1
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 (
24  const FiniteElement &el, ElementTransformation &Trans,
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,
33  ElementTransformation &Trans, DenseMatrix &elmat )
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,
49  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
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 (
66  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
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,
75  ElementTransformation &Trans, DenseMatrix &elmat)
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 (
92  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
93 {
94  bfi -> AssembleElementMatrix (el, Trans, elmat);
95  elmat.Lump();
96 }
97 
98 void InverseIntegrator::AssembleElementMatrix(
99  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
100 {
101  integrator->AssembleElementMatrix(el, Trans, elmat);
102  elmat.Invert();
103 }
104 
105 void SumIntegrator::AssembleElementMatrix(
106  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
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 ConvectionIntegrator::AssembleElementMatrix(
573  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
574 {
575  int nd = el.GetDof();
576  int dim = el.GetDim();
577 
578  elmat.SetSize(nd);
579  dshape.SetSize(nd,dim);
580  adjJ.SetSize(dim);
581  shape.SetSize(nd);
582  vec2.SetSize(dim);
583  BdFidxT.SetSize(nd);
584 
585  Vector vec1;
586 
587  const IntegrationRule *ir = IntRule;
588  if (ir == NULL)
589  {
590  int order = Trans.OrderGrad(&el) + Trans.Order() + el.GetOrder();
591  ir = &IntRules.Get(el.GetGeomType(), order);
592  }
593 
594  Q.Eval(Q_ir, Trans, *ir);
595 
596  elmat = 0.0;
597  for (int i = 0; i < ir->GetNPoints(); i++)
598  {
599  const IntegrationPoint &ip = ir->IntPoint(i);
600  el.CalcDShape(ip, dshape);
601  el.CalcShape(ip, shape);
602 
603  Trans.SetIntPoint(&ip);
604  CalcAdjugate(Trans.Jacobian(), adjJ);
605  Q_ir.GetColumnReference(i, vec1);
606  vec1 *= alpha * ip.weight;
607 
608  adjJ.Mult(vec1, vec2);
609  dshape.Mult(vec2, BdFidxT);
610 
611  AddMultVWt(shape, BdFidxT, elmat);
612  }
613 }
614 
615 
616 void GroupConvectionIntegrator::AssembleElementMatrix(
617  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
618 {
619  int nd = el.GetDof();
620  int dim = el.GetDim();
621 
622  elmat.SetSize(nd);
623  dshape.SetSize(nd,dim);
624  adjJ.SetSize(dim);
625  shape.SetSize(nd);
626  grad.SetSize(nd,dim);
627 
628  const IntegrationRule *ir = IntRule;
629  if (ir == NULL)
630  {
631  int order = Trans.OrderGrad(&el) + el.GetOrder();
632  ir = &IntRules.Get(el.GetGeomType(), order);
633  }
634 
635  Q.Eval(Q_nodal, Trans, el.GetNodes()); // sets the size of Q_nodal
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 
647  Mult(dshape, adjJ, grad);
648 
649  double w = alpha * ip.weight;
650 
651  // elmat(k,l) += \sum_s w*shape(k)*Q_nodal(s,k)*grad(l,s)
652  for (int k = 0; k < nd; k++)
653  {
654  double wsk = w*shape(k);
655  for (int l = 0; l < nd; l++)
656  {
657  double a = 0.0;
658  for (int s = 0; s < dim; s++)
659  {
660  a += Q_nodal(s,k)*grad(l,s);
661  }
662  elmat(k,l) += wsk*a;
663  }
664  }
665  }
666 }
667 
668 
669 void VectorMassIntegrator::AssembleElementMatrix
671  DenseMatrix &elmat )
672 {
673  int nd = el.GetDof();
674  int spaceDim = Trans.GetSpaceDim();
675 
676  double norm;
677 
678  // Get vdim from VQ, MQ, or the space dimension
679  int vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : spaceDim);
680 
681  elmat.SetSize(nd*vdim);
682  shape.SetSize(nd);
683  partelmat.SetSize(nd);
684  if (VQ)
685  {
686  vec.SetSize(vdim);
687  }
688  else if (MQ)
689  {
690  mcoeff.SetSize(vdim);
691  }
692 
693  const IntegrationRule *ir = IntRule;
694  if (ir == NULL)
695  {
696  int order = 2 * el.GetOrder() + Trans.OrderW() + Q_order;
697 
698  if (el.Space() == FunctionSpace::rQk)
699  {
700  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
701  }
702  else
703  {
704  ir = &IntRules.Get(el.GetGeomType(), order);
705  }
706  }
707 
708  elmat = 0.0;
709  for (int s = 0; s < ir->GetNPoints(); s++)
710  {
711  const IntegrationPoint &ip = ir->IntPoint(s);
712  el.CalcShape(ip, shape);
713 
714  Trans.SetIntPoint (&ip);
715  norm = ip.weight * Trans.Weight();
716 
717  MultVVt(shape, partelmat);
718 
719  if (VQ)
720  {
721  VQ->Eval(vec, Trans, ip);
722  for (int k = 0; k < vdim; k++)
723  {
724  elmat.AddMatrix(norm*vec(k), partelmat, nd*k, nd*k);
725  }
726  }
727  else if (MQ)
728  {
729  MQ->Eval(mcoeff, Trans, ip);
730  for (int i = 0; i < vdim; i++)
731  for (int j = 0; j < vdim; j++)
732  {
733  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, nd*i, nd*j);
734  }
735  }
736  else
737  {
738  if (Q)
739  {
740  norm *= Q->Eval(Trans, ip);
741  }
742  partelmat *= norm;
743  for (int k = 0; k < vdim; k++)
744  {
745  elmat.AddMatrix(partelmat, nd*k, nd*k);
746  }
747  }
748  }
749 }
750 
751 void VectorMassIntegrator::AssembleElementMatrix2(
752  const FiniteElement &trial_fe, const FiniteElement &test_fe,
753  ElementTransformation &Trans, DenseMatrix &elmat)
754 {
755  int tr_nd = trial_fe.GetDof();
756  int te_nd = test_fe.GetDof();
757  int dim = trial_fe.GetDim();
758  int vdim;
759 
760  double norm;
761 
762  // Get vdim from the ElementTransformation Trans ?
763  vdim = (VQ) ? (VQ -> GetVDim()) : ((MQ) ? (MQ -> GetVDim()) : (dim));
764 
765  elmat.SetSize(te_nd*vdim, tr_nd*vdim);
766  shape.SetSize(tr_nd);
767  te_shape.SetSize(te_nd);
768  partelmat.SetSize(te_nd, tr_nd);
769  if (VQ)
770  {
771  vec.SetSize(vdim);
772  }
773  else if (MQ)
774  {
775  mcoeff.SetSize(vdim);
776  }
777 
778  const IntegrationRule *ir = IntRule;
779  if (ir == NULL)
780  {
781  int order = (trial_fe.GetOrder() + test_fe.GetOrder() +
782  Trans.OrderW() + Q_order);
783 
784  if (trial_fe.Space() == FunctionSpace::rQk)
785  {
786  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
787  }
788  else
789  {
790  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
791  }
792  }
793 
794  elmat = 0.0;
795  for (int s = 0; s < ir->GetNPoints(); s++)
796  {
797  const IntegrationPoint &ip = ir->IntPoint(s);
798  trial_fe.CalcShape(ip, shape);
799  test_fe.CalcShape(ip, te_shape);
800 
801  Trans.SetIntPoint(&ip);
802  norm = ip.weight * Trans.Weight();
803 
804  MultVWt(te_shape, shape, partelmat);
805 
806  if (VQ)
807  {
808  VQ->Eval(vec, Trans, ip);
809  for (int k = 0; k < vdim; k++)
810  {
811  elmat.AddMatrix(norm*vec(k), partelmat, te_nd*k, tr_nd*k);
812  }
813  }
814  else if (MQ)
815  {
816  MQ->Eval(mcoeff, Trans, ip);
817  for (int i = 0; i < vdim; i++)
818  for (int j = 0; j < vdim; j++)
819  {
820  elmat.AddMatrix(norm*mcoeff(i,j), partelmat, te_nd*i, tr_nd*j);
821  }
822  }
823  else
824  {
825  if (Q)
826  {
827  norm *= Q->Eval(Trans, ip);
828  }
829  partelmat *= norm;
830  for (int k = 0; k < vdim; k++)
831  {
832  elmat.AddMatrix(partelmat, te_nd*k, tr_nd*k);
833  }
834  }
835  }
836 }
837 
838 void VectorFEDivergenceIntegrator::AssembleElementMatrix2(
839  const FiniteElement &trial_fe, const FiniteElement &test_fe,
840  ElementTransformation &Trans, DenseMatrix &elmat)
841 {
842  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
843 
844 #ifdef MFEM_THREAD_SAFE
845  Vector divshape(trial_nd), shape(test_nd);
846 #else
847  divshape.SetSize(trial_nd);
848  shape.SetSize(test_nd);
849 #endif
850 
851  elmat.SetSize(test_nd, trial_nd);
852 
853  const IntegrationRule *ir = IntRule;
854  if (ir == NULL)
855  {
856  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
857  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
858  }
859 
860  elmat = 0.0;
861  for (i = 0; i < ir->GetNPoints(); i++)
862  {
863  const IntegrationPoint &ip = ir->IntPoint(i);
864  trial_fe.CalcDivShape(ip, divshape);
865  test_fe.CalcShape(ip, shape);
866  double w = ip.weight;
867  if (Q)
868  {
869  Trans.SetIntPoint(&ip);
870  w *= Q->Eval(Trans, ip);
871  }
872  shape *= w;
873  AddMultVWt(shape, divshape, elmat);
874  }
875 }
876 
877 void VectorFECurlIntegrator::AssembleElementMatrix2(
878  const FiniteElement &trial_fe, const FiniteElement &test_fe,
879  ElementTransformation &Trans, DenseMatrix &elmat)
880 {
881  int trial_nd = trial_fe.GetDof(), test_nd = test_fe.GetDof(), i;
882  int dim = trial_fe.GetDim();
883 
884 #ifdef MFEM_THREAD_SAFE
885  DenseMatrix curlshapeTrial(trial_nd, dim);
886  DenseMatrix curlshapeTrial_dFT(trial_nd, dim);
887  DenseMatrix vshapeTest(test_nd, dim);
888 #else
889  curlshapeTrial.SetSize(trial_nd, dim);
890  curlshapeTrial_dFT.SetSize(trial_nd, dim);
891  vshapeTest.SetSize(test_nd, dim);
892 #endif
893 
894  elmat.SetSize(test_nd, trial_nd);
895 
896  const IntegrationRule *ir = IntRule;
897  if (ir == NULL)
898  {
899  int order = trial_fe.GetOrder() + test_fe.GetOrder() - 1; // <--
900  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
901  }
902 
903  elmat = 0.0;
904  for (i = 0; i < ir->GetNPoints(); i++)
905  {
906  const IntegrationPoint &ip = ir->IntPoint(i);
907 
908  Trans.SetIntPoint(&ip);
909  trial_fe.CalcCurlShape(ip, curlshapeTrial);
910  MultABt(curlshapeTrial, Trans.Jacobian(), curlshapeTrial_dFT);
911  test_fe.CalcVShape(Trans, vshapeTest);
912  double w = ip.weight;
913 
914  if (Q)
915  {
916  w *= Q->Eval(Trans, ip);
917  }
918  vshapeTest *= w;
919  AddMultABt(vshapeTest, curlshapeTrial_dFT, elmat);
920  }
921 }
922 
923 void DerivativeIntegrator::AssembleElementMatrix2 (
924  const FiniteElement &trial_fe,
925  const FiniteElement &test_fe,
926  ElementTransformation &Trans,
927  DenseMatrix &elmat)
928 {
929  int dim = trial_fe.GetDim();
930  int trial_nd = trial_fe.GetDof();
931  int test_nd = test_fe.GetDof();
932 
933  int i, l;
934  double det;
935 
936  elmat.SetSize (test_nd,trial_nd);
937  dshape.SetSize (trial_nd,dim);
938  dshapedxt.SetSize(trial_nd,dim);
939  dshapedxi.SetSize(trial_nd);
940  invdfdx.SetSize(dim);
941  shape.SetSize (test_nd);
942 
943  const IntegrationRule *ir = IntRule;
944  if (ir == NULL)
945  {
946  int order;
947  if (trial_fe.Space() == FunctionSpace::Pk)
948  {
949  order = trial_fe.GetOrder() + test_fe.GetOrder() - 1;
950  }
951  else
952  {
953  order = trial_fe.GetOrder() + test_fe.GetOrder() + dim;
954  }
955 
956  if (trial_fe.Space() == FunctionSpace::rQk)
957  {
958  ir = &RefinedIntRules.Get(trial_fe.GetGeomType(), order);
959  }
960  else
961  {
962  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
963  }
964  }
965 
966  elmat = 0.0;
967  for (i = 0; i < ir->GetNPoints(); i++)
968  {
969  const IntegrationPoint &ip = ir->IntPoint(i);
970 
971  trial_fe.CalcDShape(ip, dshape);
972 
973  Trans.SetIntPoint (&ip);
974  CalcInverse (Trans.Jacobian(), invdfdx);
975  det = Trans.Weight();
976  Mult (dshape, invdfdx, dshapedxt);
977 
978  test_fe.CalcShape(ip, shape);
979 
980  for (l = 0; l < trial_nd; l++)
981  {
982  dshapedxi(l) = dshapedxt(l,xi);
983  }
984 
985  shape *= Q.Eval(Trans,ip) * det * ip.weight;
986  AddMultVWt (shape, dshapedxi, elmat);
987  }
988 }
989 
990 void CurlCurlIntegrator::AssembleElementMatrix
992  DenseMatrix &elmat )
993 {
994  int nd = el.GetDof();
995  int dim = el.GetDim();
996  int dimc = (dim == 3) ? 3 : 1;
997  double w;
998 
999 #ifdef MFEM_THREAD_SAFE
1000  DenseMatrix curlshape(nd,dimc), curlshape_dFt(nd,dimc);
1001 #else
1002  curlshape.SetSize(nd,dimc);
1003  curlshape_dFt.SetSize(nd,dimc);
1004 #endif
1005  elmat.SetSize(nd);
1006 
1007  const IntegrationRule *ir = IntRule;
1008  if (ir == NULL)
1009  {
1010  int order;
1011  if (el.Space() == FunctionSpace::Pk)
1012  {
1013  order = 2*el.GetOrder() - 2;
1014  }
1015  else
1016  {
1017  order = 2*el.GetOrder();
1018  }
1019 
1020  ir = &IntRules.Get(el.GetGeomType(), order);
1021  }
1022 
1023  elmat = 0.0;
1024  for (int i = 0; i < ir->GetNPoints(); i++)
1025  {
1026  const IntegrationPoint &ip = ir->IntPoint(i);
1027 
1028  Trans.SetIntPoint (&ip);
1029 
1030  w = ip.weight / Trans.Weight();
1031 
1032  if ( dim == 3 )
1033  {
1034  el.CalcCurlShape(ip, curlshape);
1035  MultABt(curlshape, Trans.Jacobian(), curlshape_dFt);
1036  }
1037  else
1038  {
1039  el.CalcCurlShape(ip, curlshape_dFt);
1040  }
1041 
1042  if (Q)
1043  {
1044  w *= Q->Eval(Trans, ip);
1045  }
1046 
1047  AddMult_a_AAt(w, curlshape_dFt, elmat);
1048  }
1049 }
1050 
1051 void CurlCurlIntegrator
1052 ::ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans,
1053  Vector &u, const FiniteElement &fluxelem, Vector &flux,
1054  int with_coef)
1055 {
1056 #ifdef MFEM_THREAD_SAFE
1057  DenseMatrix projcurl;
1058 #endif
1059 
1060  fluxelem.ProjectCurl(el, Trans, projcurl);
1061 
1062  flux.SetSize(projcurl.Height());
1063  projcurl.Mult(u, flux);
1064 
1065  // TODO: Q, wcoef?
1066 }
1067 
1068 double CurlCurlIntegrator::ComputeFluxEnergy(const FiniteElement &fluxelem,
1069  ElementTransformation &Trans,
1070  Vector &flux, Vector *d_energy)
1071 {
1072  int nd = fluxelem.GetDof();
1073  int dim = fluxelem.GetDim();
1074 
1075 #ifdef MFEM_THREAD_SAFE
1076  DenseMatrix vshape;
1077 #endif
1078  vshape.SetSize(nd, dim);
1079  pointflux.SetSize(dim);
1080  if (d_energy) { vec.SetSize(dim); }
1081 
1082  int order = 2 * fluxelem.GetOrder(); // <--
1083  const IntegrationRule &ir = IntRules.Get(fluxelem.GetGeomType(), order);
1084 
1085  double energy = 0.0;
1086  if (d_energy) { *d_energy = 0.0; }
1087 
1088  Vector* pfluxes = NULL;
1089  if (d_energy)
1090  {
1091  pfluxes = new Vector[ir.GetNPoints()];
1092  }
1093 
1094  for (int i = 0; i < ir.GetNPoints(); i++)
1095  {
1096  const IntegrationPoint &ip = ir.IntPoint(i);
1097  Trans.SetIntPoint(&ip);
1098 
1099  fluxelem.CalcVShape(Trans, vshape);
1100  // fluxelem.CalcVShape(ip, vshape);
1101  vshape.MultTranspose(flux, pointflux);
1102 
1103  double w = Trans.Weight() * ip.weight;
1104 
1105  double e = w * (pointflux * pointflux);
1106 
1107  if (Q)
1108  {
1109  // TODO
1110  }
1111 
1112  energy += e;
1113 
1114 #if ANISO_EXPERIMENTAL
1115  if (d_energy)
1116  {
1117  pfluxes[i].SetSize(dim);
1118  Trans.Jacobian().MultTranspose(pointflux, pfluxes[i]);
1119 
1120  /*
1121  DenseMatrix Jadj(dim, dim);
1122  CalcAdjugate(Trans.Jacobian(), Jadj);
1123  pfluxes[i].SetSize(dim);
1124  Jadj.Mult(pointflux, pfluxes[i]);
1125  */
1126 
1127  // pfluxes[i] = pointflux;
1128  }
1129 #endif
1130  }
1131 
1132  if (d_energy)
1133  {
1134 #if ANISO_EXPERIMENTAL
1135  *d_energy = 0.0;
1136  Vector tmp;
1137 
1138  int n = (int) round(pow(ir.GetNPoints(), 1.0/3.0));
1139  MFEM_ASSERT(n*n*n == ir.GetNPoints(), "");
1140 
1141  // hack: get total variation of 'pointflux' in the x,y,z directions
1142  for (int k = 0; k < n; k++)
1143  for (int l = 0; l < n; l++)
1144  for (int m = 0; m < n; m++)
1145  {
1146  Vector &vec = pfluxes[(k*n + l)*n + m];
1147  if (m > 0)
1148  {
1149  tmp = vec; tmp -= pfluxes[(k*n + l)*n + (m-1)];
1150  (*d_energy)[0] += (tmp * tmp);
1151  }
1152  if (l > 0)
1153  {
1154  tmp = vec; tmp -= pfluxes[(k*n + (l-1))*n + m];
1155  (*d_energy)[1] += (tmp * tmp);
1156  }
1157  if (k > 0)
1158  {
1159  tmp = vec; tmp -= pfluxes[((k-1)*n + l)*n + m];
1160  (*d_energy)[2] += (tmp * tmp);
1161  }
1162  }
1163 #else
1164  *d_energy = 1.0;
1165 #endif
1166 
1167  delete [] pfluxes;
1168  }
1169 
1170  return energy;
1171 }
1172 
1173 void VectorCurlCurlIntegrator::AssembleElementMatrix(
1174  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1175 {
1176  int dim = el.GetDim();
1177  int dof = el.GetDof();
1178  int cld = (dim*(dim-1))/2;
1179 
1180 #ifdef MFEM_THREAD_SAFE
1181  DenseMatrix dshape_hat(dof, dim), dshape(dof, dim);
1182  DenseMatrix curlshape(dim*dof, cld), Jadj(dim);
1183 #else
1184  dshape_hat.SetSize(dof, dim);
1185  dshape.SetSize(dof, dim);
1186  curlshape.SetSize(dim*dof, cld);
1187  Jadj.SetSize(dim);
1188 #endif
1189 
1190  const IntegrationRule *ir = IntRule;
1191  if (ir == NULL)
1192  {
1193  // use the same integration rule as diffusion
1194  int order = 2 * Trans.OrderGrad(&el);
1195  ir = &IntRules.Get(el.GetGeomType(), order);
1196  }
1197 
1198  elmat = 0.0;
1199  for (int i = 0; i < ir->GetNPoints(); i++)
1200  {
1201  const IntegrationPoint &ip = ir->IntPoint(i);
1202  el.CalcDShape(ip, dshape_hat);
1203 
1204  Trans.SetIntPoint(&ip);
1205  CalcAdjugate(Trans.Jacobian(), Jadj);
1206  double w = ip.weight / Trans.Weight();
1207 
1208  Mult(dshape_hat, Jadj, dshape);
1209  dshape.GradToCurl(curlshape);
1210 
1211  if (Q)
1212  {
1213  w *= Q->Eval(Trans, ip);
1214  }
1215 
1216  AddMult_a_AAt(w, curlshape, elmat);
1217  }
1218 }
1219 
1220 double VectorCurlCurlIntegrator::GetElementEnergy(
1221  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
1222 {
1223  int dim = el.GetDim();
1224  int dof = el.GetDof();
1225 
1226 #ifdef MFEM_THREAD_SAFE
1227  DenseMatrix dshape_hat(dof, dim), Jadj(dim), grad_hat(dim), grad(dim);
1228 #else
1229  dshape_hat.SetSize(dof, dim);
1230 
1231  Jadj.SetSize(dim);
1232  grad_hat.SetSize(dim);
1233  grad.SetSize(dim);
1234 #endif
1235  DenseMatrix elfun_mat(elfun.GetData(), dof, dim);
1236 
1237  const IntegrationRule *ir = IntRule;
1238  if (ir == NULL)
1239  {
1240  // use the same integration rule as diffusion
1241  int order = 2 * Tr.OrderGrad(&el);
1242  ir = &IntRules.Get(el.GetGeomType(), order);
1243  }
1244 
1245  double energy = 0.;
1246  for (int i = 0; i < ir->GetNPoints(); i++)
1247  {
1248  const IntegrationPoint &ip = ir->IntPoint(i);
1249  el.CalcDShape(ip, dshape_hat);
1250 
1251  MultAtB(elfun_mat, dshape_hat, grad_hat);
1252 
1253  Tr.SetIntPoint(&ip);
1254  CalcAdjugate(Tr.Jacobian(), Jadj);
1255  double w = ip.weight / Tr.Weight();
1256 
1257  Mult(grad_hat, Jadj, grad);
1258 
1259  if (dim == 2)
1260  {
1261  double curl = grad(0,1) - grad(1,0);
1262  w *= curl * curl;
1263  }
1264  else
1265  {
1266  double curl_x = grad(2,1) - grad(1,2);
1267  double curl_y = grad(0,2) - grad(2,0);
1268  double curl_z = grad(1,0) - grad(0,1);
1269  w *= curl_x * curl_x + curl_y * curl_y + curl_z * curl_z;
1270  }
1271 
1272  if (Q)
1273  {
1274  w *= Q->Eval(Tr, ip);
1275  }
1276 
1277  energy += w;
1278  }
1279 
1280  elfun_mat.ClearExternalData();
1281 
1282  return 0.5 * energy;
1283 }
1284 
1285 
1286 void VectorFEMassIntegrator::AssembleElementMatrix(
1287  const FiniteElement &el,
1288  ElementTransformation &Trans,
1289  DenseMatrix &elmat)
1290 {
1291  int dof = el.GetDof();
1292  int spaceDim = Trans.GetSpaceDim();
1293 
1294  double w;
1295 
1296 #ifdef MFEM_THREAD_SAFE
1297  Vector D(VQ ? VQ->GetVDim() : 0);
1298  DenseMatrix trial_vshape(dof, spaceDim);
1299  DenseMatrix K(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1300 #else
1301  trial_vshape.SetSize(dof, spaceDim);
1302  D.SetSize(VQ ? VQ->GetVDim() : 0);
1303  K.SetSize(MQ ? MQ->GetVDim() : 0, MQ ? MQ->GetVDim() : 0);
1304 #endif
1305  DenseMatrix tmp(trial_vshape.Height(), K.Width());
1306 
1307  elmat.SetSize(dof);
1308  elmat = 0.0;
1309 
1310  const IntegrationRule *ir = IntRule;
1311  if (ir == NULL)
1312  {
1313  // int order = 2 * el.GetOrder();
1314  int order = Trans.OrderW() + 2 * el.GetOrder();
1315  ir = &IntRules.Get(el.GetGeomType(), order);
1316  }
1317 
1318  for (int i = 0; i < ir->GetNPoints(); i++)
1319  {
1320  const IntegrationPoint &ip = ir->IntPoint(i);
1321 
1322  Trans.SetIntPoint (&ip);
1323 
1324  el.CalcVShape(Trans, trial_vshape);
1325 
1326  w = ip.weight * Trans.Weight();
1327  if (MQ)
1328  {
1329  MQ->Eval(K, Trans, ip);
1330  K *= w;
1331  Mult(trial_vshape,K,tmp);
1332  AddMultABt(tmp,trial_vshape,elmat);
1333  }
1334  else if (VQ)
1335  {
1336  VQ->Eval(D, Trans, ip);
1337  D *= w;
1338  AddMultADAt(trial_vshape, D, elmat);
1339  }
1340  else
1341  {
1342  if (Q)
1343  {
1344  w *= Q -> Eval (Trans, ip);
1345  }
1346  AddMult_a_AAt (w, trial_vshape, elmat);
1347  }
1348  }
1349 }
1350 
1351 void VectorFEMassIntegrator::AssembleElementMatrix2(
1352  const FiniteElement &trial_fe, const FiniteElement &test_fe,
1353  ElementTransformation &Trans, DenseMatrix &elmat)
1354 {
1355  if ( test_fe.GetRangeType() == FiniteElement::SCALAR && VQ )
1356  {
1357  // assume test_fe is scalar FE and trial_fe is vector FE
1358  int dim = test_fe.GetDim();
1359  int trial_dof = trial_fe.GetDof();
1360  int test_dof = test_fe.GetDof();
1361  double w;
1362 
1363  if (MQ)
1364  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1365  " is not implemented for tensor materials");
1366 
1367 #ifdef MFEM_THREAD_SAFE
1368  DenseMatrix trial_vshape(trial_dof, dim);
1369  Vector shape(test_dof);
1370  Vector D(dim);
1371 #else
1372  trial_vshape.SetSize(trial_dof, dim);
1373  shape.SetSize(test_dof);
1374  D.SetSize(dim);
1375 #endif
1376 
1377  elmat.SetSize (test_dof, trial_dof);
1378 
1379  const IntegrationRule *ir = IntRule;
1380  if (ir == NULL)
1381  {
1382  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1383  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1384  }
1385 
1386  elmat = 0.0;
1387  for (int i = 0; i < ir->GetNPoints(); i++)
1388  {
1389  const IntegrationPoint &ip = ir->IntPoint(i);
1390 
1391  Trans.SetIntPoint (&ip);
1392 
1393  trial_fe.CalcVShape(Trans, trial_vshape);
1394  test_fe.CalcShape(ip, shape);
1395 
1396  w = ip.weight * Trans.Weight();
1397  VQ->Eval(D, Trans, ip);
1398  D *= w;
1399 
1400  for (int d = 0; d < dim; d++)
1401  {
1402  for (int j = 0; j < test_dof; j++)
1403  {
1404  for (int k = 0; k < trial_dof; k++)
1405  {
1406  elmat(j, k) += D[d] * shape(j) * trial_vshape(k, d);
1407  }
1408  }
1409  }
1410  }
1411  }
1412  else if ( test_fe.GetRangeType() == FiniteElement::SCALAR )
1413  {
1414  // assume test_fe is scalar FE and trial_fe is vector FE
1415  int dim = test_fe.GetDim();
1416  int trial_dof = trial_fe.GetDof();
1417  int test_dof = test_fe.GetDof();
1418  double w;
1419 
1420  if (VQ || MQ)
1421  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1422  " is not implemented for vector/tensor permeability");
1423 
1424 #ifdef MFEM_THREAD_SAFE
1425  DenseMatrix trial_vshape(trial_dof, dim);
1426  Vector shape(test_dof);
1427 #else
1428  trial_vshape.SetSize(trial_dof, dim);
1429  shape.SetSize(test_dof);
1430 #endif
1431 
1432  elmat.SetSize (dim*test_dof, trial_dof);
1433 
1434  const IntegrationRule *ir = IntRule;
1435  if (ir == NULL)
1436  {
1437  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1438  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1439  }
1440 
1441  elmat = 0.0;
1442  for (int i = 0; i < ir->GetNPoints(); i++)
1443  {
1444  const IntegrationPoint &ip = ir->IntPoint(i);
1445 
1446  Trans.SetIntPoint (&ip);
1447 
1448  trial_fe.CalcVShape(Trans, trial_vshape);
1449  test_fe.CalcShape(ip, shape);
1450 
1451  w = ip.weight * Trans.Weight();
1452  if (Q)
1453  {
1454  w *= Q -> Eval (Trans, ip);
1455  }
1456 
1457  for (int d = 0; d < dim; d++)
1458  {
1459  for (int j = 0; j < test_dof; j++)
1460  {
1461  for (int k = 0; k < trial_dof; k++)
1462  {
1463  elmat(d * test_dof + j, k) += w * shape(j) * trial_vshape(k, d);
1464  }
1465  }
1466  }
1467  }
1468  }
1469  else
1470  {
1471  // assume both test_fe and trial_fe are vector FE
1472  int dim = test_fe.GetDim();
1473  int trial_dof = trial_fe.GetDof();
1474  int test_dof = test_fe.GetDof();
1475  double w;
1476 
1477  if (VQ || MQ)
1478  mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
1479  " is not implemented for vector/tensor permeability");
1480 
1481 #ifdef MFEM_THREAD_SAFE
1482  DenseMatrix trial_vshape(trial_dof, dim);
1483  DenseMatrix test_vshape(test_dof,dim);
1484 #else
1485  trial_vshape.SetSize(trial_dof, dim);
1486  test_vshape.SetSize(test_dof,dim);
1487 #endif
1488 
1489  elmat.SetSize (test_dof, trial_dof);
1490 
1491  const IntegrationRule *ir = IntRule;
1492  if (ir == NULL)
1493  {
1494  int order = (Trans.OrderW() + test_fe.GetOrder() + trial_fe.GetOrder());
1495  ir = &IntRules.Get(test_fe.GetGeomType(), order);
1496  }
1497 
1498  elmat = 0.0;
1499  for (int i = 0; i < ir->GetNPoints(); i++)
1500  {
1501  const IntegrationPoint &ip = ir->IntPoint(i);
1502 
1503  Trans.SetIntPoint (&ip);
1504 
1505  trial_fe.CalcVShape(Trans, trial_vshape);
1506  test_fe.CalcVShape(Trans, test_vshape);
1507 
1508  w = ip.weight * Trans.Weight();
1509  if (Q)
1510  {
1511  w *= Q -> Eval (Trans, ip);
1512  }
1513 
1514  for (int d = 0; d < dim; d++)
1515  {
1516  for (int j = 0; j < test_dof; j++)
1517  {
1518  for (int k = 0; k < trial_dof; k++)
1519  {
1520  elmat(j, k) += w * test_vshape(j, d) * trial_vshape(k, d);
1521  }
1522  }
1523  }
1524  }
1525  }
1526 }
1527 
1528 void VectorDivergenceIntegrator::AssembleElementMatrix2(
1529  const FiniteElement &trial_fe,
1530  const FiniteElement &test_fe,
1531  ElementTransformation &Trans,
1532  DenseMatrix &elmat)
1533 {
1534  int dim = trial_fe.GetDim();
1535  int trial_dof = trial_fe.GetDof();
1536  int test_dof = test_fe.GetDof();
1537  double c;
1538 
1539  dshape.SetSize (trial_dof, dim);
1540  gshape.SetSize (trial_dof, dim);
1541  Jadj.SetSize (dim);
1542  divshape.SetSize (dim*trial_dof);
1543  shape.SetSize (test_dof);
1544 
1545  elmat.SetSize (test_dof, dim*trial_dof);
1546 
1547  const IntegrationRule *ir = IntRule;
1548  if (ir == NULL)
1549  {
1550  int order = Trans.OrderGrad(&trial_fe) + test_fe.GetOrder();
1551  ir = &IntRules.Get(trial_fe.GetGeomType(), order);
1552  }
1553 
1554  elmat = 0.0;
1555 
1556  for (int i = 0; i < ir -> GetNPoints(); i++)
1557  {
1558  const IntegrationPoint &ip = ir->IntPoint(i);
1559 
1560  trial_fe.CalcDShape (ip, dshape);
1561  test_fe.CalcShape (ip, shape);
1562 
1563  Trans.SetIntPoint (&ip);
1564  CalcAdjugate(Trans.Jacobian(), Jadj);
1565 
1566  Mult (dshape, Jadj, gshape);
1567 
1568  gshape.GradToDiv (divshape);
1569 
1570  c = ip.weight;
1571  if (Q)
1572  {
1573  c *= Q -> Eval (Trans, ip);
1574  }
1575 
1576  // elmat += c * shape * divshape ^ t
1577  shape *= c;
1578  AddMultVWt (shape, divshape, elmat);
1579  }
1580 }
1581 
1582 
1583 void DivDivIntegrator::AssembleElementMatrix(
1584  const FiniteElement &el,
1585  ElementTransformation &Trans,
1586  DenseMatrix &elmat)
1587 {
1588  int dof = el.GetDof();
1589  double c;
1590 
1591 #ifdef MFEM_THREAD_SAFE
1592  Vector divshape(dof);
1593 #else
1594  divshape.SetSize(dof);
1595 #endif
1596  elmat.SetSize(dof);
1597 
1598  const IntegrationRule *ir = IntRule;
1599  if (ir == NULL)
1600  {
1601  int order = 2 * el.GetOrder() - 2; // <--- OK for RTk
1602  ir = &IntRules.Get(el.GetGeomType(), order);
1603  }
1604 
1605  elmat = 0.0;
1606 
1607  for (int i = 0; i < ir -> GetNPoints(); i++)
1608  {
1609  const IntegrationPoint &ip = ir->IntPoint(i);
1610 
1611  el.CalcDivShape (ip, divshape);
1612 
1613  Trans.SetIntPoint (&ip);
1614  c = ip.weight / Trans.Weight();
1615 
1616  if (Q)
1617  {
1618  c *= Q -> Eval (Trans, ip);
1619  }
1620 
1621  // elmat += c * divshape * divshape ^ t
1622  AddMult_a_VVt (c, divshape, elmat);
1623  }
1624 }
1625 
1626 
1627 void VectorDiffusionIntegrator::AssembleElementMatrix(
1628  const FiniteElement &el,
1629  ElementTransformation &Trans,
1630  DenseMatrix &elmat)
1631 {
1632  int dim = el.GetDim();
1633  int dof = el.GetDof();
1634 
1635  double norm;
1636 
1637  elmat.SetSize (dim * dof);
1638 
1639  Jinv. SetSize (dim);
1640  dshape.SetSize (dof, dim);
1641  gshape.SetSize (dof, dim);
1642  pelmat.SetSize (dof);
1643 
1644  const IntegrationRule *ir = IntRule;
1645  if (ir == NULL)
1646  {
1647  // integrand is rational function if det(J) is not constant
1648  int order = 2 * Trans.OrderGrad(&el); // order of the numerator
1649  if (el.Space() == FunctionSpace::rQk)
1650  {
1651  ir = &RefinedIntRules.Get(el.GetGeomType(), order);
1652  }
1653  else
1654  {
1655  ir = &IntRules.Get(el.GetGeomType(), order);
1656  }
1657  }
1658 
1659  elmat = 0.0;
1660 
1661  for (int i = 0; i < ir -> GetNPoints(); i++)
1662  {
1663  const IntegrationPoint &ip = ir->IntPoint(i);
1664 
1665  el.CalcDShape (ip, dshape);
1666 
1667  Trans.SetIntPoint (&ip);
1668  norm = ip.weight * Trans.Weight();
1669  CalcInverse (Trans.Jacobian(), Jinv);
1670 
1671  Mult (dshape, Jinv, gshape);
1672 
1673  MultAAt (gshape, pelmat);
1674 
1675  if (Q)
1676  {
1677  norm *= Q -> Eval (Trans, ip);
1678  }
1679 
1680  pelmat *= norm;
1681 
1682  for (int d = 0; d < dim; d++)
1683  {
1684  for (int k = 0; k < dof; k++)
1685  for (int l = 0; l < dof; l++)
1686  {
1687  elmat (dof*d+k, dof*d+l) += pelmat (k, l);
1688  }
1689  }
1690  }
1691 }
1692 
1693 
1694 void ElasticityIntegrator::AssembleElementMatrix(
1695  const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
1696 {
1697  int dof = el.GetDof();
1698  int dim = el.GetDim();
1699  double w, L, M;
1700 
1701 #ifdef MFEM_THREAD_SAFE
1702  DenseMatrix dshape(dof, dim), Jinv(dim), gshape(dof, dim), pelmat(dof);
1703  Vector divshape(dim*dof);
1704 #else
1705  Jinv.SetSize(dim);
1706  dshape.SetSize(dof, dim);
1707  gshape.SetSize(dof, dim);
1708  pelmat.SetSize(dof);
1709  divshape.SetSize(dim*dof);
1710 #endif
1711 
1712  elmat.SetSize(dof * dim);
1713 
1714  const IntegrationRule *ir = IntRule;
1715  if (ir == NULL)
1716  {
1717  int order = 2 * Trans.OrderGrad(&el); // correct order?
1718  ir = &IntRules.Get(el.GetGeomType(), order);
1719  }
1720 
1721  elmat = 0.0;
1722 
1723  for (int i = 0; i < ir -> GetNPoints(); i++)
1724  {
1725  const IntegrationPoint &ip = ir->IntPoint(i);
1726 
1727  el.CalcDShape(ip, dshape);
1728 
1729  Trans.SetIntPoint(&ip);
1730  w = ip.weight * Trans.Weight();
1731  CalcInverse(Trans.Jacobian(), Jinv);
1732  Mult(dshape, Jinv, gshape);
1733  MultAAt(gshape, pelmat);
1734  gshape.GradToDiv (divshape);
1735 
1736  M = mu->Eval(Trans, ip);
1737  if (lambda)
1738  {
1739  L = lambda->Eval(Trans, ip);
1740  }
1741  else
1742  {
1743  L = q_lambda * M;
1744  M = q_mu * M;
1745  }
1746 
1747  if (L != 0.0)
1748  {
1749  AddMult_a_VVt(L * w, divshape, elmat);
1750  }
1751 
1752  if (M != 0.0)
1753  {
1754  for (int d = 0; d < dim; d++)
1755  {
1756  for (int k = 0; k < dof; k++)
1757  for (int l = 0; l < dof; l++)
1758  {
1759  elmat (dof*d+k, dof*d+l) += (M * w) * pelmat(k, l);
1760  }
1761  }
1762  for (int i = 0; i < dim; i++)
1763  for (int j = 0; j < dim; j++)
1764  {
1765  for (int k = 0; k < dof; k++)
1766  for (int l = 0; l < dof; l++)
1767  elmat(dof*i+k, dof*j+l) +=
1768  (M * w) * gshape(k, j) * gshape(l, i);
1769  // + (L * w) * gshape(k, i) * gshape(l, j)
1770  }
1771  }
1772  }
1773 }
1774 
1775 void DGTraceIntegrator::AssembleFaceMatrix(const FiniteElement &el1,
1776  const FiniteElement &el2,
1778  DenseMatrix &elmat)
1779 {
1780  int dim, ndof1, ndof2;
1781 
1782  double un, a, b, w;
1783 
1784  dim = el1.GetDim();
1785  ndof1 = el1.GetDof();
1786  Vector vu(dim), nor(dim);
1787 
1788  if (Trans.Elem2No >= 0)
1789  {
1790  ndof2 = el2.GetDof();
1791  }
1792  else
1793  {
1794  ndof2 = 0;
1795  }
1796 
1797  shape1.SetSize(ndof1);
1798  shape2.SetSize(ndof2);
1799  elmat.SetSize(ndof1 + ndof2);
1800  elmat = 0.0;
1801 
1802  const IntegrationRule *ir = IntRule;
1803  if (ir == NULL)
1804  {
1805  int order;
1806  // Assuming order(u)==order(mesh)
1807  if (Trans.Elem2No >= 0)
1808  order = (min(Trans.Elem1->OrderW(), Trans.Elem2->OrderW()) +
1809  2*max(el1.GetOrder(), el2.GetOrder()));
1810  else
1811  {
1812  order = Trans.Elem1->OrderW() + 2*el1.GetOrder();
1813  }
1814  if (el1.Space() == FunctionSpace::Pk)
1815  {
1816  order++;
1817  }
1818  ir = &IntRules.Get(Trans.FaceGeom, order);
1819  }
1820 
1821  for (int p = 0; p < ir->GetNPoints(); p++)
1822  {
1823  const IntegrationPoint &ip = ir->IntPoint(p);
1824  IntegrationPoint eip1, eip2;
1825  Trans.Loc1.Transform(ip, eip1);
1826  if (ndof2)
1827  {
1828  Trans.Loc2.Transform(ip, eip2);
1829  }
1830  el1.CalcShape(eip1, shape1);
1831 
1832  Trans.Face->SetIntPoint(&ip);
1833  Trans.Elem1->SetIntPoint(&eip1);
1834 
1835  u->Eval(vu, *Trans.Elem1, eip1);
1836 
1837  if (dim == 1)
1838  {
1839  nor(0) = 2*eip1.x - 1.0;
1840  }
1841  else
1842  {
1843  CalcOrtho(Trans.Face->Jacobian(), nor);
1844  }
1845 
1846  un = vu * nor;
1847  a = 0.5 * alpha * un;
1848  b = beta * fabs(un);
1849  // note: if |alpha/2|==|beta| then |a|==|b|, i.e. (a==b) or (a==-b)
1850  // and therefore two blocks in the element matrix contribution
1851  // (from the current quadrature point) are 0
1852 
1853  if (rho)
1854  {
1855  double rho_p;
1856  if (un >= 0.0 && ndof2)
1857  {
1858  Trans.Elem2->SetIntPoint(&eip2);
1859  rho_p = rho->Eval(*Trans.Elem2, eip2);
1860  }
1861  else
1862  {
1863  rho_p = rho->Eval(*Trans.Elem1, eip1);
1864  }
1865  a *= rho_p;
1866  b *= rho_p;
1867  }
1868 
1869  w = ip.weight * (a+b);
1870  if (w != 0.0)
1871  {
1872  for (int i = 0; i < ndof1; i++)
1873  for (int j = 0; j < ndof1; j++)
1874  {
1875  elmat(i, j) += w * shape1(i) * shape1(j);
1876  }
1877  }
1878 
1879  if (ndof2)
1880  {
1881  el2.CalcShape(eip2, shape2);
1882 
1883  if (w != 0.0)
1884  for (int i = 0; i < ndof2; i++)
1885  for (int j = 0; j < ndof1; j++)
1886  {
1887  elmat(ndof1+i, j) -= w * shape2(i) * shape1(j);
1888  }
1889 
1890  w = ip.weight * (b-a);
1891  if (w != 0.0)
1892  {
1893  for (int i = 0; i < ndof2; i++)
1894  for (int j = 0; j < ndof2; j++)
1895  {
1896  elmat(ndof1+i, ndof1+j) += w * shape2(i) * shape2(j);
1897  }
1898 
1899  for (int i = 0; i < ndof1; i++)
1900  for (int j = 0; j < ndof2; j++)
1901  {
1902  elmat(i, ndof1+j) -= w * shape1(i) * shape2(j);
1903  }
1904  }
1905  }
1906  }
1907 }
1908 
1909 void DGDiffusionIntegrator::AssembleFaceMatrix(
1910  const FiniteElement &el1, const FiniteElement &el2,
1911  FaceElementTransformations &Trans, DenseMatrix &elmat)
1912 {
1913  int dim, ndof1, ndof2, ndofs;
1914  bool kappa_is_nonzero = (kappa != 0.);
1915  double w, wq = 0.0;
1916 
1917  dim = el1.GetDim();
1918  ndof1 = el1.GetDof();
1919 
1920  nor.SetSize(dim);
1921  nh.SetSize(dim);
1922  ni.SetSize(dim);
1923  adjJ.SetSize(dim);
1924  if (MQ)
1925  {
1926  mq.SetSize(dim);
1927  }
1928 
1929  shape1.SetSize(ndof1);
1930  dshape1.SetSize(ndof1, dim);
1931  dshape1dn.SetSize(ndof1);
1932  if (Trans.Elem2No >= 0)
1933  {
1934  ndof2 = el2.GetDof();
1935  shape2.SetSize(ndof2);
1936  dshape2.SetSize(ndof2, dim);
1937  dshape2dn.SetSize(ndof2);
1938  }
1939  else
1940  {
1941  ndof2 = 0;
1942  }
1943 
1944  ndofs = ndof1 + ndof2;
1945  elmat.SetSize(ndofs);
1946  elmat = 0.0;
1947  if (kappa_is_nonzero)
1948  {
1949  jmat.SetSize(ndofs);
1950  jmat = 0.;
1951  }
1952 
1953  const IntegrationRule *ir = IntRule;
1954  if (ir == NULL)
1955  {
1956  // a simple choice for the integration order; is this OK?
1957  int order;
1958  if (ndof2)
1959  {
1960  order = 2*max(el1.GetOrder(), el2.GetOrder());
1961  }
1962  else
1963  {
1964  order = 2*el1.GetOrder();
1965  }
1966  ir = &IntRules.Get(Trans.FaceGeom, order);
1967  }
1968 
1969  // assemble: < {(Q \nabla u).n},[v] > --> elmat
1970  // kappa < {h^{-1} Q} [u],[v] > --> jmat
1971  for (int p = 0; p < ir->GetNPoints(); p++)
1972  {
1973  const IntegrationPoint &ip = ir->IntPoint(p);
1974  IntegrationPoint eip1, eip2;
1975 
1976  Trans.Loc1.Transform(ip, eip1);
1977  Trans.Face->SetIntPoint(&ip);
1978  if (dim == 1)
1979  {
1980  nor(0) = 2*eip1.x - 1.0;
1981  }
1982  else
1983  {
1984  CalcOrtho(Trans.Face->Jacobian(), nor);
1985  }
1986 
1987  el1.CalcShape(eip1, shape1);
1988  el1.CalcDShape(eip1, dshape1);
1989  Trans.Elem1->SetIntPoint(&eip1);
1990  w = ip.weight/Trans.Elem1->Weight();
1991  if (ndof2)
1992  {
1993  w /= 2;
1994  }
1995  if (!MQ)
1996  {
1997  if (Q)
1998  {
1999  w *= Q->Eval(*Trans.Elem1, eip1);
2000  }
2001  ni.Set(w, nor);
2002  }
2003  else
2004  {
2005  nh.Set(w, nor);
2006  MQ->Eval(mq, *Trans.Elem1, eip1);
2007  mq.MultTranspose(nh, ni);
2008  }
2009  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
2010  adjJ.Mult(ni, nh);
2011  if (kappa_is_nonzero)
2012  {
2013  wq = ni * nor;
2014  }
2015  // Note: in the jump term, we use 1/h1 = |nor|/det(J1) which is
2016  // independent of Loc1 and always gives the size of element 1 in
2017  // direction perpendicular to the face. Indeed, for linear transformation
2018  // |nor|=measure(face)/measure(ref. face),
2019  // det(J1)=measure(element)/measure(ref. element),
2020  // and the ratios measure(ref. element)/measure(ref. face) are
2021  // compatible for all element/face pairs.
2022  // For example: meas(ref. tetrahedron)/meas(ref. triangle) = 1/3, and
2023  // for any tetrahedron vol(tet)=(1/3)*height*area(base).
2024  // For interior faces: q_e/h_e=(q1/h1+q2/h2)/2.
2025 
2026  dshape1.Mult(nh, dshape1dn);
2027  for (int i = 0; i < ndof1; i++)
2028  for (int j = 0; j < ndof1; j++)
2029  {
2030  elmat(i, j) += shape1(i) * dshape1dn(j);
2031  }
2032 
2033  if (ndof2)
2034  {
2035  Trans.Loc2.Transform(ip, eip2);
2036  el2.CalcShape(eip2, shape2);
2037  el2.CalcDShape(eip2, dshape2);
2038  Trans.Elem2->SetIntPoint(&eip2);
2039  w = ip.weight/2/Trans.Elem2->Weight();
2040  if (!MQ)
2041  {
2042  if (Q)
2043  {
2044  w *= Q->Eval(*Trans.Elem2, eip2);
2045  }
2046  ni.Set(w, nor);
2047  }
2048  else
2049  {
2050  nh.Set(w, nor);
2051  MQ->Eval(mq, *Trans.Elem2, eip2);
2052  mq.MultTranspose(nh, ni);
2053  }
2054  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
2055  adjJ.Mult(ni, nh);
2056  if (kappa_is_nonzero)
2057  {
2058  wq += ni * nor;
2059  }
2060 
2061  dshape2.Mult(nh, dshape2dn);
2062 
2063  for (int i = 0; i < ndof1; i++)
2064  for (int j = 0; j < ndof2; j++)
2065  {
2066  elmat(i, ndof1 + j) += shape1(i) * dshape2dn(j);
2067  }
2068 
2069  for (int i = 0; i < ndof2; i++)
2070  for (int j = 0; j < ndof1; j++)
2071  {
2072  elmat(ndof1 + i, j) -= shape2(i) * dshape1dn(j);
2073  }
2074 
2075  for (int i = 0; i < ndof2; i++)
2076  for (int j = 0; j < ndof2; j++)
2077  {
2078  elmat(ndof1 + i, ndof1 + j) -= shape2(i) * dshape2dn(j);
2079  }
2080  }
2081 
2082  if (kappa_is_nonzero)
2083  {
2084  // only assemble the lower triangular part of jmat
2085  wq *= kappa;
2086  for (int i = 0; i < ndof1; i++)
2087  {
2088  const double wsi = wq*shape1(i);
2089  for (int j = 0; j <= i; j++)
2090  {
2091  jmat(i, j) += wsi * shape1(j);
2092  }
2093  }
2094  if (ndof2)
2095  {
2096  for (int i = 0; i < ndof2; i++)
2097  {
2098  const int i2 = ndof1 + i;
2099  const double wsi = wq*shape2(i);
2100  for (int j = 0; j < ndof1; j++)
2101  {
2102  jmat(i2, j) -= wsi * shape1(j);
2103  }
2104  for (int j = 0; j <= i; j++)
2105  {
2106  jmat(i2, ndof1 + j) += wsi * shape2(j);
2107  }
2108  }
2109  }
2110  }
2111  }
2112 
2113  // elmat := -elmat + sigma*elmat^t + jmat
2114  if (kappa_is_nonzero)
2115  {
2116  for (int i = 0; i < ndofs; i++)
2117  {
2118  for (int j = 0; j < i; j++)
2119  {
2120  double aij = elmat(i,j), aji = elmat(j,i), mij = jmat(i,j);
2121  elmat(i,j) = sigma*aji - aij + mij;
2122  elmat(j,i) = sigma*aij - aji + mij;
2123  }
2124  elmat(i,i) = (sigma - 1.)*elmat(i,i) + jmat(i,i);
2125  }
2126  }
2127  else
2128  {
2129  for (int i = 0; i < ndofs; i++)
2130  {
2131  for (int j = 0; j < i; j++)
2132  {
2133  double aij = elmat(i,j), aji = elmat(j,i);
2134  elmat(i,j) = sigma*aji - aij;
2135  elmat(j,i) = sigma*aij - aji;
2136  }
2137  elmat(i,i) *= (sigma - 1.);
2138  }
2139  }
2140 }
2141 
2142 void TraceJumpIntegrator::AssembleFaceMatrix(
2143  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2144  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2145  DenseMatrix &elmat)
2146 {
2147  int i, j, face_ndof, ndof1, ndof2;
2148  int order;
2149 
2150  double w;
2151 
2152  face_ndof = trial_face_fe.GetDof();
2153  ndof1 = test_fe1.GetDof();
2154 
2155  face_shape.SetSize(face_ndof);
2156  shape1.SetSize(ndof1);
2157 
2158  if (Trans.Elem2No >= 0)
2159  {
2160  ndof2 = test_fe2.GetDof();
2161  shape2.SetSize(ndof2);
2162  }
2163  else
2164  {
2165  ndof2 = 0;
2166  }
2167 
2168  elmat.SetSize(ndof1 + ndof2, face_ndof);
2169  elmat = 0.0;
2170 
2171  const IntegrationRule *ir = IntRule;
2172  if (ir == NULL)
2173  {
2174  if (Trans.Elem2No >= 0)
2175  {
2176  order = max(test_fe1.GetOrder(), test_fe2.GetOrder());
2177  }
2178  else
2179  {
2180  order = test_fe1.GetOrder();
2181  }
2182  order += trial_face_fe.GetOrder();
2183  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2184  {
2185  order += Trans.Face->OrderW();
2186  }
2187  ir = &IntRules.Get(Trans.FaceGeom, order);
2188  }
2189 
2190  for (int p = 0; p < ir->GetNPoints(); p++)
2191  {
2192  const IntegrationPoint &ip = ir->IntPoint(p);
2193  IntegrationPoint eip1, eip2;
2194  // Trace finite element shape function
2195  Trans.Face->SetIntPoint(&ip);
2196  trial_face_fe.CalcShape(ip, face_shape);
2197  // Side 1 finite element shape function
2198  Trans.Loc1.Transform(ip, eip1);
2199  test_fe1.CalcShape(eip1, shape1);
2200  Trans.Elem1->SetIntPoint(&eip1);
2201  if (ndof2)
2202  {
2203  // Side 2 finite element shape function
2204  Trans.Loc2.Transform(ip, eip2);
2205  test_fe2.CalcShape(eip2, shape2);
2206  Trans.Elem2->SetIntPoint(&eip2);
2207  }
2208  w = ip.weight;
2209  if (trial_face_fe.GetMapType() == FiniteElement::VALUE)
2210  {
2211  w *= Trans.Face->Weight();
2212  }
2213  face_shape *= w;
2214  for (i = 0; i < ndof1; i++)
2215  for (j = 0; j < face_ndof; j++)
2216  {
2217  elmat(i, j) += shape1(i) * face_shape(j);
2218  }
2219  if (ndof2)
2220  {
2221  // Subtract contribution from side 2
2222  for (i = 0; i < ndof2; i++)
2223  for (j = 0; j < face_ndof; j++)
2224  {
2225  elmat(ndof1+i, j) -= shape2(i) * face_shape(j);
2226  }
2227  }
2228  }
2229 }
2230 
2231 void NormalTraceJumpIntegrator::AssembleFaceMatrix(
2232  const FiniteElement &trial_face_fe, const FiniteElement &test_fe1,
2233  const FiniteElement &test_fe2, FaceElementTransformations &Trans,
2234  DenseMatrix &elmat)
2235 {
2236  int i, j, face_ndof, ndof1, ndof2, dim;
2237  int order;
2238 
2239  MFEM_VERIFY(trial_face_fe.GetMapType() == FiniteElement::VALUE, "");
2240 
2241  face_ndof = trial_face_fe.GetDof();
2242  ndof1 = test_fe1.GetDof();
2243  dim = test_fe1.GetDim();
2244 
2245  face_shape.SetSize(face_ndof);
2246  normal.SetSize(dim);
2247  shape1.SetSize(ndof1,dim);
2248  shape1_n.SetSize(ndof1);
2249 
2250  if (Trans.Elem2No >= 0)
2251  {
2252  ndof2 = test_fe2.GetDof();
2253  shape2.SetSize(ndof2,dim);
2254  shape2_n.SetSize(ndof2);
2255  }
2256  else
2257  {
2258  ndof2 = 0;
2259  }
2260 
2261  elmat.SetSize(ndof1 + ndof2, face_ndof);
2262  elmat = 0.0;
2263 
2264  const IntegrationRule *ir = IntRule;
2265  if (ir == NULL)
2266  {
2267  if (Trans.Elem2No >= 0)
2268  {
2269  order = max(test_fe1.GetOrder(), test_fe2.GetOrder()) - 1;
2270  }
2271  else
2272  {
2273  order = test_fe1.GetOrder() - 1;
2274  }
2275  order += trial_face_fe.GetOrder();
2276  ir = &IntRules.Get(Trans.FaceGeom, order);
2277  }
2278 
2279  for (int p = 0; p < ir->GetNPoints(); p++)
2280  {
2281  const IntegrationPoint &ip = ir->IntPoint(p);
2282  IntegrationPoint eip1, eip2;
2283  // Trace finite element shape function
2284  trial_face_fe.CalcShape(ip, face_shape);
2285  Trans.Loc1.Transf.SetIntPoint(&ip);
2286  CalcOrtho(Trans.Loc1.Transf.Jacobian(), normal);
2287  // Side 1 finite element shape function
2288  Trans.Loc1.Transform(ip, eip1);
2289  test_fe1.CalcVShape(eip1, shape1);
2290  shape1.Mult(normal, shape1_n);
2291  if (ndof2)
2292  {
2293  // Side 2 finite element shape function
2294  Trans.Loc2.Transform(ip, eip2);
2295  test_fe2.CalcVShape(eip2, shape2);
2296  Trans.Loc2.Transf.SetIntPoint(&ip);
2297  CalcOrtho(Trans.Loc2.Transf.Jacobian(), normal);
2298  shape2.Mult(normal, shape2_n);
2299  }
2300  face_shape *= ip.weight;
2301  for (i = 0; i < ndof1; i++)
2302  for (j = 0; j < face_ndof; j++)
2303  {
2304  elmat(i, j) -= shape1_n(i) * face_shape(j);
2305  }
2306  if (ndof2)
2307  {
2308  // Subtract contribution from side 2
2309  for (i = 0; i < ndof2; i++)
2310  for (j = 0; j < face_ndof; j++)
2311  {
2312  elmat(ndof1+i, j) += shape2_n(i) * face_shape(j);
2313  }
2314  }
2315  }
2316 }
2317 
2318 
2319 void NormalInterpolator::AssembleElementMatrix2(
2320  const FiniteElement &dom_fe, const FiniteElement &ran_fe,
2321  ElementTransformation &Trans, DenseMatrix &elmat)
2322 {
2323  int spaceDim = Trans.GetSpaceDim();
2324  elmat.SetSize(ran_fe.GetDof(), spaceDim*dom_fe.GetDof());
2325  Vector n(spaceDim), shape(dom_fe.GetDof());
2326 
2327  const IntegrationRule &ran_nodes = ran_fe.GetNodes();
2328  for (int i = 0; i < ran_nodes.Size(); i++)
2329  {
2330  const IntegrationPoint &ip = ran_nodes.IntPoint(i);
2331  Trans.SetIntPoint(&ip);
2332  CalcOrtho(Trans.Jacobian(), n);
2333  dom_fe.CalcShape(ip, shape);
2334  for (int j = 0; j < shape.Size(); j++)
2335  {
2336  for (int d = 0; d < spaceDim; d++)
2337  {
2338  elmat(i, j+d*shape.Size()) = shape(j)*n(d);
2339  }
2340  }
2341  }
2342 }
2343 
2344 }
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:3182
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:82
ElementTransformation * Face
Definition: eltrans.hpp:119
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3489
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:3467
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:235
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:259
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:445
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:33
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2864
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:120
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:119
double * GetData() const
Definition: vector.hpp:88
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3100
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:48
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:120
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:568
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:3532
virtual double Weight()=0
virtual const DenseMatrix & Jacobian()
Definition: eltrans.cpp:52
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2972
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:85
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:2284
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:3458
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:3305
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:3124
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:2535
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:110
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:2348
ElementTransformation * Elem1
Definition: eltrans.hpp:119
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:3362
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:3138
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:3419