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