MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lininteg.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 
13 #include <cmath>
14 #include "fem.hpp"
15 
16 namespace mfem
17 {
18 
20  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
21 {
22  mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
23 }
24 
25 
28  Vector &elvect)
29 {
30  int dof = el.GetDof();
31 
32  shape.SetSize(dof); // vector of size dof
33  elvect.SetSize(dof);
34  elvect = 0.0;
35 
36  const IntegrationRule *ir = IntRule;
37  if (ir == NULL)
38  {
39  // ir = &IntRules.Get(el.GetGeomType(),
40  // oa * el.GetOrder() + ob + Tr.OrderW());
41  ir = &IntRules.Get(el.GetGeomType(), oa * el.GetOrder() + ob);
42  }
43 
44  for (int i = 0; i < ir->GetNPoints(); i++)
45  {
46  const IntegrationPoint &ip = ir->IntPoint(i);
47 
48  Tr.SetIntPoint (&ip);
49  double val = Tr.Weight() * Q.Eval(Tr, ip);
50 
51  el.CalcShape(ip, shape);
52 
53  add(elvect, ip.weight * val, shape, elvect);
54  }
55 }
56 
58  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
59 {
60  int dof = el.GetDof();
61 
62  shape.SetSize(dof); // vector of size dof
63  elvect.SetSize(dof);
64  elvect = 0.0;
65 
66  const IntegrationRule *ir = IntRule;
67  if (ir == NULL)
68  {
69  int intorder = oa * el.GetOrder() + ob; // <----------
70  ir = &IntRules.Get(el.GetGeomType(), intorder);
71  }
72 
73  for (int i = 0; i < ir->GetNPoints(); i++)
74  {
75  const IntegrationPoint &ip = ir->IntPoint(i);
76 
77  Tr.SetIntPoint (&ip);
78  double val = Tr.Weight() * Q.Eval(Tr, ip);
79 
80  el.CalcShape(ip, shape);
81 
82  add(elvect, ip.weight * val, shape, elvect);
83  }
84 }
85 
87  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
88 {
89  int dim = el.GetDim()+1;
90  int dof = el.GetDof();
91  Vector nor(dim), Qvec;
92 
93  shape.SetSize(dof);
94  elvect.SetSize(dof);
95  elvect = 0.0;
96 
97  const IntegrationRule *ir = IntRule;
98  if (ir == NULL)
99  {
100  int intorder = oa * el.GetOrder() + ob; // <----------
101  ir = &IntRules.Get(el.GetGeomType(), intorder);
102  }
103 
104  for (int i = 0; i < ir->GetNPoints(); i++)
105  {
106  const IntegrationPoint &ip = ir->IntPoint(i);
107 
108  Tr.SetIntPoint(&ip);
109  CalcOrtho(Tr.Jacobian(), nor);
110  Q.Eval(Qvec, Tr, ip);
111 
112  el.CalcShape(ip, shape);
113 
114  elvect.Add(ip.weight*(Qvec*nor), shape);
115  }
116 }
117 
119  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
120 {
121  int dim = el.GetDim()+1;
122  int dof = el.GetDof();
123  Vector tangent(dim), Qvec;
124 
125  shape.SetSize(dof);
126  elvect.SetSize(dof);
127  elvect = 0.0;
128 
129  if (dim != 2)
130  {
131  mfem_error("These methods make sense only in 2D problems.");
132  }
133 
134  const IntegrationRule *ir = IntRule;
135  if (ir == NULL)
136  {
137  int intorder = oa * el.GetOrder() + ob; // <----------
138  ir = &IntRules.Get(el.GetGeomType(), intorder);
139  }
140 
141  for (int i = 0; i < ir->GetNPoints(); i++)
142  {
143  const IntegrationPoint &ip = ir->IntPoint(i);
144 
145  Tr.SetIntPoint(&ip);
146  const DenseMatrix &Jac = Tr.Jacobian();
147  tangent(0) = Jac(0,0);
148  tangent(1) = Jac(1,0);
149 
150  Q.Eval(Qvec, Tr, ip);
151 
152  el.CalcShape(ip, shape);
153 
154  add(elvect, ip.weight*(Qvec*tangent), shape, elvect);
155  }
156 }
157 
159  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
160 {
161  int vdim = Q.GetVDim();
162  int dof = el.GetDof();
163 
164  double val,cf;
165 
166  shape.SetSize(dof); // vector of size dof
167 
168  elvect.SetSize(dof * vdim);
169  elvect = 0.0;
170 
171  const IntegrationRule *ir = IntRule;
172  if (ir == NULL)
173  {
174  int intorder = el.GetOrder() + 1;
175  ir = &IntRules.Get(el.GetGeomType(), intorder);
176  }
177 
178  for (int i = 0; i < ir->GetNPoints(); i++)
179  {
180  const IntegrationPoint &ip = ir->IntPoint(i);
181 
182  Tr.SetIntPoint (&ip);
183  val = Tr.Weight();
184 
185  el.CalcShape(ip, shape);
186  Q.Eval (Qvec, Tr, ip);
187 
188  for (int k = 0; k < vdim; k++)
189  {
190  cf = val * Qvec(k);
191 
192  for (int s = 0; s < dof; s++)
193  {
194  elvect(dof*k+s) += ip.weight * cf * shape(s);
195  }
196  }
197  }
198 }
199 
201  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
202 {
203  int vdim = Q.GetVDim();
204  int dof = el.GetDof();
205 
206  shape.SetSize(dof);
207  vec.SetSize(vdim);
208 
209  elvect.SetSize(dof * vdim);
210  elvect = 0.0;
211 
212  const IntegrationRule *ir = IntRule;
213  if (ir == NULL)
214  {
215  int intorder = el.GetOrder() + 1;
216  ir = &IntRules.Get(el.GetGeomType(), intorder);
217  }
218 
219  for (int i = 0; i < ir->GetNPoints(); i++)
220  {
221  const IntegrationPoint &ip = ir->IntPoint(i);
222 
223  Q.Eval(vec, Tr, ip);
224  Tr.SetIntPoint (&ip);
225  vec *= Tr.Weight() * ip.weight;
226  el.CalcShape(ip, shape);
227  for (int k = 0; k < vdim; k++)
228  for (int s = 0; s < dof; s++)
229  {
230  elvect(dof*k+s) += vec(k) * shape(s);
231  }
232  }
233 }
234 
236  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
237 {
238  int vdim = Q.GetVDim();
239  int dof = el.GetDof();
240 
241  shape.SetSize(dof);
242  vec.SetSize(vdim);
243 
244  elvect.SetSize(dof * vdim);
245  elvect = 0.0;
246 
247  const IntegrationRule *ir = IntRule;
248  if (ir == NULL)
249  {
250  int intorder = el.GetOrder() + 1;
251  ir = &IntRules.Get(Tr.FaceGeom, intorder);
252  }
253 
254  for (int i = 0; i < ir->GetNPoints(); i++)
255  {
256  const IntegrationPoint &ip = ir->IntPoint(i);
257  IntegrationPoint eip;
258  Tr.Loc1.Transform(ip, eip);
259 
260  Tr.Face->SetIntPoint(&ip);
261  Q.Eval(vec, *Tr.Face, ip);
262  vec *= Tr.Face->Weight() * ip.weight;
263  el.CalcShape(eip, shape);
264  for (int k = 0; k < vdim; k++)
265  {
266  for (int s = 0; s < dof; s++)
267  {
268  elvect(dof*k+s) += vec(k) * shape(s);
269  }
270  }
271  }
272 }
273 
274 
276  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
277 {
278  int dof = el.GetDof();
279  int spaceDim = Tr.GetSpaceDim();
280 
281  vshape.SetSize(dof,spaceDim);
282  vec.SetSize(spaceDim);
283 
284  elvect.SetSize(dof);
285  elvect = 0.0;
286 
287  const IntegrationRule *ir = IntRule;
288  if (ir == NULL)
289  {
290  // int intorder = 2*el.GetOrder() - 1; // ok for O(h^{k+1}) conv. in L2
291  int intorder = 2*el.GetOrder();
292  ir = &IntRules.Get(el.GetGeomType(), intorder);
293  }
294 
295  for (int i = 0; i < ir->GetNPoints(); i++)
296  {
297  const IntegrationPoint &ip = ir->IntPoint(i);
298 
299  Tr.SetIntPoint (&ip);
300  el.CalcVShape(Tr, vshape);
301 
302  QF.Eval (vec, Tr, ip);
303  vec *= ip.weight * Tr.Weight();
304 
305  vshape.AddMult (vec, elvect);
306  }
307 }
308 
309 
311  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
312 {
313  int dim = el.GetDim()+1;
314  int dof = el.GetDof();
315 
316  shape.SetSize (dof);
317  nor.SetSize (dim);
318  elvect.SetSize (dim*dof);
319 
320  const IntegrationRule *ir = IntRule;
321  if (ir == NULL)
322  {
323  ir = &IntRules.Get(el.GetGeomType(), el.GetOrder() + 1);
324  }
325 
326  elvect = 0.0;
327  for (int i = 0; i < ir->GetNPoints(); i++)
328  {
329  const IntegrationPoint &ip = ir->IntPoint(i);
330  Tr.SetIntPoint (&ip);
331  CalcOrtho(Tr.Jacobian(), nor);
332  el.CalcShape (ip, shape);
333  nor *= Sign * ip.weight * F -> Eval (Tr, ip);
334  for (int j = 0; j < dof; j++)
335  for (int k = 0; k < dim; k++)
336  {
337  elvect(dof*k+j) += nor(k) * shape(j);
338  }
339  }
340 }
341 
342 
344  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
345 {
346  int dof = el.GetDof();
347 
348  shape.SetSize(dof);
349  elvect.SetSize(dof);
350  elvect = 0.0;
351 
352  const IntegrationRule *ir = IntRule;
353  if (ir == NULL)
354  {
355  int intorder = 2*el.GetOrder(); // <----------
356  ir = &IntRules.Get(el.GetGeomType(), intorder);
357  }
358 
359  for (int i = 0; i < ir->GetNPoints(); i++)
360  {
361  const IntegrationPoint &ip = ir->IntPoint(i);
362 
363  Tr.SetIntPoint (&ip);
364  double val = ip.weight*F.Eval(Tr, ip);
365 
366  el.CalcShape(ip, shape);
367 
368  add(elvect, val, shape, elvect);
369  }
370 }
371 
372 
374  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
375 {
376  int dof = el.GetDof();
377  DenseMatrix vshape(dof, 2);
378  Vector f_loc(3);
379  Vector f_hat(2);
380 
381  elvect.SetSize(dof);
382  elvect = 0.0;
383 
384  const IntegrationRule *ir = IntRule;
385  if (ir == NULL)
386  {
387  int intorder = 2*el.GetOrder(); // <----------
388  ir = &IntRules.Get(el.GetGeomType(), intorder);
389  }
390 
391  for (int i = 0; i < ir->GetNPoints(); i++)
392  {
393  const IntegrationPoint &ip = ir->IntPoint(i);
394 
395  Tr.SetIntPoint(&ip);
396  f.Eval(f_loc, Tr, ip);
397  Tr.Jacobian().MultTranspose(f_loc, f_hat);
398  el.CalcVShape(ip, vshape);
399 
400  Swap<double>(f_hat(0), f_hat(1));
401  f_hat(0) = -f_hat(0);
402  f_hat *= ip.weight;
403  vshape.AddMult(f_hat, elvect);
404  }
405 }
406 
407 
409  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
410 {
411  mfem_error("BoundaryFlowIntegrator::AssembleRHSElementVect\n"
412  " is not implemented as boundary integrator!\n"
413  " Use LinearForm::AddBdrFaceIntegrator instead of\n"
414  " LinearForm::AddBoundaryIntegrator.");
415 }
416 
418  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
419 {
420  int dim, ndof, order;
421  double un, w, vu_data[3], nor_data[3];
422 
423  dim = el.GetDim();
424  ndof = el.GetDof();
425  Vector vu(vu_data, dim), nor(nor_data, dim);
426 
427  const IntegrationRule *ir = IntRule;
428  if (ir == NULL)
429  {
430  // Assuming order(u)==order(mesh)
431  order = Tr.Elem1->OrderW() + 2*el.GetOrder();
432  if (el.Space() == FunctionSpace::Pk)
433  {
434  order++;
435  }
436  ir = &IntRules.Get(Tr.FaceGeom, order);
437  }
438 
439  shape.SetSize(ndof);
440  elvect.SetSize(ndof);
441  elvect = 0.0;
442 
443  for (int p = 0; p < ir->GetNPoints(); p++)
444  {
445  const IntegrationPoint &ip = ir->IntPoint(p);
446  IntegrationPoint eip;
447  Tr.Loc1.Transform(ip, eip);
448  el.CalcShape(eip, shape);
449 
450  Tr.Face->SetIntPoint(&ip);
451 
452  u->Eval(vu, *Tr.Elem1, eip);
453 
454  if (dim == 1)
455  {
456  nor(0) = 2*eip.x - 1.0;
457  }
458  else
459  {
460  CalcOrtho(Tr.Face->Jacobian(), nor);
461  }
462 
463  un = vu * nor;
464  w = 0.5*alpha*un - beta*fabs(un);
465  w *= ip.weight*f->Eval(*Tr.Elem1, eip);
466  elvect.Add(w, shape);
467  }
468 }
469 
470 
472  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
473 {
474  mfem_error("DGDirichletLFIntegrator::AssembleRHSElementVect");
475 }
476 
478  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
479 {
480  int dim, ndof;
481  bool kappa_is_nonzero = (kappa != 0.);
482  double w;
483 
484  dim = el.GetDim();
485  ndof = el.GetDof();
486 
487  nor.SetSize(dim);
488  nh.SetSize(dim);
489  ni.SetSize(dim);
490  adjJ.SetSize(dim);
491  if (MQ)
492  {
493  mq.SetSize(dim);
494  }
495 
496  shape.SetSize(ndof);
497  dshape.SetSize(ndof, dim);
498  dshape_dn.SetSize(ndof);
499 
500  elvect.SetSize(ndof);
501  elvect = 0.0;
502 
503  const IntegrationRule *ir = IntRule;
504  if (ir == NULL)
505  {
506  // a simple choice for the integration order; is this OK?
507  int order = 2*el.GetOrder();
508  ir = &IntRules.Get(Tr.FaceGeom, order);
509  }
510 
511  for (int p = 0; p < ir->GetNPoints(); p++)
512  {
513  const IntegrationPoint &ip = ir->IntPoint(p);
514  IntegrationPoint eip;
515 
516  Tr.Loc1.Transform(ip, eip);
517  Tr.Face->SetIntPoint(&ip);
518  if (dim == 1)
519  {
520  nor(0) = 2*eip.x - 1.0;
521  }
522  else
523  {
524  CalcOrtho(Tr.Face->Jacobian(), nor);
525  }
526 
527  el.CalcShape(eip, shape);
528  el.CalcDShape(eip, dshape);
529  Tr.Elem1->SetIntPoint(&eip);
530  // compute uD through the face transformation
531  w = ip.weight * uD->Eval(*Tr.Face, ip) / Tr.Elem1->Weight();
532  if (!MQ)
533  {
534  if (Q)
535  {
536  w *= Q->Eval(*Tr.Elem1, eip);
537  }
538  ni.Set(w, nor);
539  }
540  else
541  {
542  nh.Set(w, nor);
543  MQ->Eval(mq, *Tr.Elem1, eip);
544  mq.MultTranspose(nh, ni);
545  }
547  adjJ.Mult(ni, nh);
548 
550  elvect.Add(sigma, dshape_dn);
551 
552  if (kappa_is_nonzero)
553  {
554  elvect.Add(kappa*(ni*nor), shape);
555  }
556  }
557 }
558 
559 
561  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
562 {
563  mfem_error("DGElasticityDirichletLFIntegrator::AssembleRHSElementVect");
564 }
565 
567  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
568 {
569  MFEM_ASSERT(Tr.Elem2No < 0, "interior boundary is not supported");
570 
571 #ifdef MFEM_THREAD_SAFE
572  Vector shape;
576  Vector nor;
579  Vector u_dir;
580 #endif
581 
582  const int dim = el.GetDim();
583  const int ndofs = el.GetDof();
584  const int nvdofs = dim*ndofs;
585 
586  elvect.SetSize(nvdofs);
587  elvect = 0.0;
588 
589  adjJ.SetSize(dim);
590  shape.SetSize(ndofs);
591  dshape.SetSize(ndofs, dim);
592  dshape_ps.SetSize(ndofs, dim);
593  nor.SetSize(dim);
594  dshape_dn.SetSize(ndofs);
595  dshape_du.SetSize(ndofs);
596  u_dir.SetSize(dim);
597 
598  const IntegrationRule *ir = IntRule;
599  if (ir == NULL)
600  {
601  const int order = 2*el.GetOrder(); // <-----
602  ir = &IntRules.Get(Tr.FaceGeom, order);
603  }
604 
605  for (int pi = 0; pi < ir->GetNPoints(); ++pi)
606  {
607  const IntegrationPoint &ip = ir->IntPoint(pi);
608  IntegrationPoint eip;
609  Tr.Loc1.Transform(ip, eip);
610  Tr.Face->SetIntPoint(&ip);
611  Tr.Elem1->SetIntPoint(&eip);
612 
613  // Evaluate the Dirichlet b.c. using the face transformation.
614  uD.Eval(u_dir, *Tr.Face, ip);
615 
616  el.CalcShape(eip, shape);
617  el.CalcDShape(eip, dshape);
618 
620  Mult(dshape, adjJ, dshape_ps);
621 
622  if (dim == 1)
623  {
624  nor(0) = 2*eip.x - 1.0;
625  }
626  else
627  {
628  CalcOrtho(Tr.Face->Jacobian(), nor);
629  }
630 
631  double wL, wM, jcoef;
632  {
633  const double w = ip.weight / Tr.Elem1->Weight();
634  wL = w * lambda->Eval(*Tr.Elem1, eip);
635  wM = w * mu->Eval(*Tr.Elem1, eip);
636  jcoef = kappa * (wL + 2.0*wM) * (nor*nor);
637  dshape_ps.Mult(nor, dshape_dn);
638  dshape_ps.Mult(u_dir, dshape_du);
639  }
640 
641  // alpha < uD, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
642  // + kappa < h^{-1} (lambda + 2 mu) uD, v >
643 
644  // i = idof + ndofs * im
645  // v_phi(i,d) = delta(im,d) phi(idof)
646  // div(v_phi(i)) = dphi(idof,im)
647  // (grad(v_phi(i)))(k,l) = delta(im,k) dphi(idof,l)
648  //
649  // term 1:
650  // alpha < uD, lambda div(v_phi(i)) n >
651  // alpha lambda div(v_phi(i)) (uD.n) =
652  // alpha lambda dphi(idof,im) (uD.n) --> quadrature -->
653  // ip.weight/det(J1) alpha lambda (uD.nor) dshape_ps(idof,im) =
654  // alpha * wL * (u_dir*nor) * dshape_ps(idof,im)
655  // term 2:
656  // < alpha uD, mu grad(v_phi(i)).n > =
657  // alpha mu uD^T grad(v_phi(i)) n =
658  // alpha mu uD(k) delta(im,k) dphi(idof,l) n(l) =
659  // alpha mu uD(im) dphi(idof,l) n(l) --> quadrature -->
660  // ip.weight/det(J1) alpha mu uD(im) dshape_ps(idof,l) nor(l) =
661  // alpha * wM * u_dir(im) * dshape_dn(idof)
662  // term 3:
663  // < alpha uD, mu (grad(v_phi(i)))^T n > =
664  // alpha mu n^T grad(v_phi(i)) uD =
665  // alpha mu n(k) delta(im,k) dphi(idof,l) uD(l) =
666  // alpha mu n(im) dphi(idof,l) uD(l) --> quadrature -->
667  // ip.weight/det(J1) alpha mu nor(im) dshape_ps(idof,l) uD(l) =
668  // alpha * wM * nor(im) * dshape_du(idof)
669  // term j:
670  // < kappa h^{-1} (lambda + 2 mu) uD, v_phi(i) > =
671  // kappa/h (lambda + 2 mu) uD(k) v_phi(i,k) =
672  // kappa/h (lambda + 2 mu) uD(k) delta(im,k) phi(idof) =
673  // kappa/h (lambda + 2 mu) uD(im) phi(idof) --> quadrature -->
674  // [ 1/h = |nor|/det(J1) ]
675  // ip.weight/det(J1) |nor|^2 kappa (lambda + 2 mu) uD(im) phi(idof) =
676  // jcoef * u_dir(im) * shape(idof)
677 
678  wM *= alpha;
679  const double t1 = alpha * wL * (u_dir*nor);
680  for (int im = 0, i = 0; im < dim; ++im)
681  {
682  const double t2 = wM * u_dir(im);
683  const double t3 = wM * nor(im);
684  const double tj = jcoef * u_dir(im);
685  for (int idof = 0; idof < ndofs; ++idof, ++i)
686  {
687  elvect(i) += (t1*dshape_ps(idof,im) + t2*dshape_dn(idof) +
688  t3*dshape_du(idof) + tj*shape(idof));
689  }
690  }
691  }
692 }
693 
694 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:222
Abstract class for Finite Elements.
Definition: fe.hpp:46
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:343
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:115
ElementTransformation * Face
Definition: eltrans.hpp:158
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
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:38
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:861
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:310
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:408
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:560
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:51
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2902
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:124
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:158
Polynomials of order k.
Definition: fe.hpp:34
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:127
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3132
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:57
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:264
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:225
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:193
int dim
Definition: ex3.cpp:47
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:159
const DenseMatrix & Jacobian()
Definition: eltrans.hpp:64
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:200
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...
MatrixCoefficient * MQ
Definition: lininteg.hpp:292
int GetVDim()
Returns dimension of the vector.
void AddMult(const Vector &x, Vector &y) const
y += A.x
Definition: densemat.cpp:216
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:118
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:471
int GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:118
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:121
void mfem_error(const char *msg)
Definition: error.cpp:106
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:256
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:310
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:223
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:373
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:205
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:86
ElementTransformation * Elem1
Definition: eltrans.hpp:158
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Vector data type.
Definition: vector.hpp:36
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:275
virtual int GetSpaceDim()=0
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:26
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:82
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...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:343