MFEM  v4.4.0
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-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 
13 #include "fem.hpp"
14 #include <cmath>
15 
16 namespace mfem
17 {
18 
20  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
21 {
22  mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
23 }
24 
26  const FiniteElement &el1, const FiniteElement &el2,
28 {
29  mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
30 }
31 
34  Vector &elvect)
35 {
36  int dof = el.GetDof();
37 
38  shape.SetSize(dof); // vector of size dof
39  elvect.SetSize(dof);
40  elvect = 0.0;
41 
42  const IntegrationRule *ir = IntRule;
43  if (ir == NULL)
44  {
45  // ir = &IntRules.Get(el.GetGeomType(),
46  // oa * el.GetOrder() + ob + Tr.OrderW());
47  ir = &IntRules.Get(el.GetGeomType(), oa * el.GetOrder() + ob);
48  }
49 
50  for (int i = 0; i < ir->GetNPoints(); i++)
51  {
52  const IntegrationPoint &ip = ir->IntPoint(i);
53 
54  Tr.SetIntPoint (&ip);
55  double val = Tr.Weight() * Q.Eval(Tr, ip);
56 
57  el.CalcPhysShape(Tr, shape);
58 
59  add(elvect, ip.weight * val, shape, elvect);
60  }
61 }
62 
64  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
65 {
66  MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
67  elvect.SetSize(fe.GetDof());
68  fe.CalcPhysShape(Trans, elvect);
69  elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
70 }
71 
73  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
74 {
75  int dof = el.GetDof();
76  int spaceDim = Tr.GetSpaceDim();
77 
78  dshape.SetSize(dof, spaceDim);
79 
80  elvect.SetSize(dof);
81  elvect = 0.0;
82 
83  const IntegrationRule *ir = IntRule;
84  if (ir == NULL)
85  {
86  int intorder = 2 * el.GetOrder();
87  ir = &IntRules.Get(el.GetGeomType(), intorder);
88  }
89 
90  for (int i = 0; i < ir->GetNPoints(); i++)
91  {
92  const IntegrationPoint &ip = ir->IntPoint(i);
93 
94  Tr.SetIntPoint(&ip);
95  el.CalcPhysDShape(Tr, dshape);
96 
97  Q.Eval(Qvec, Tr, ip);
98  Qvec *= ip.weight * Tr.Weight();
99 
100  dshape.AddMult(Qvec, elvect);
101  }
102 }
103 
105  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
106 {
107  MFEM_ASSERT(vec_delta != NULL,"coefficient must be VectorDeltaCoefficient");
108  int dof = fe.GetDof();
109  int spaceDim = Trans.GetSpaceDim();
110 
111  dshape.SetSize(dof, spaceDim);
112  fe.CalcPhysDShape(Trans, dshape);
113 
114  vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
115 
116  elvect.SetSize(dof);
117  dshape.Mult(Qvec, elvect);
118 }
119 
121  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
122 {
123  int dof = el.GetDof();
124 
125  shape.SetSize(dof); // vector of size dof
126  elvect.SetSize(dof);
127  elvect = 0.0;
128 
129  const IntegrationRule *ir = IntRule;
130  if (ir == NULL)
131  {
132  int intorder = oa * el.GetOrder() + ob; // <----------
133  ir = &IntRules.Get(el.GetGeomType(), intorder);
134  }
135 
136  for (int i = 0; i < ir->GetNPoints(); i++)
137  {
138  const IntegrationPoint &ip = ir->IntPoint(i);
139 
140  Tr.SetIntPoint (&ip);
141  double val = Tr.Weight() * Q.Eval(Tr, ip);
142 
143  el.CalcShape(ip, shape);
144 
145  add(elvect, ip.weight * val, shape, elvect);
146  }
147 }
148 
150  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
151 {
152  int dof = el.GetDof();
153 
154  shape.SetSize(dof); // vector of size dof
155  elvect.SetSize(dof);
156  elvect = 0.0;
157 
158  const IntegrationRule *ir = IntRule;
159  if (ir == NULL)
160  {
161  int intorder = oa * el.GetOrder() + ob; // <------ user control
162  ir = &IntRules.Get(Tr.FaceGeom, intorder); // of integration order
163  }
164 
165  for (int i = 0; i < ir->GetNPoints(); i++)
166  {
167  const IntegrationPoint &ip = ir->IntPoint(i);
168 
169  // Set the integration point in the face and the neighboring element
170  Tr.SetAllIntPoints(&ip);
171 
172  // Access the neighboring element's integration point
173  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
174 
175  double val = Tr.Face->Weight() * ip.weight * Q.Eval(*Tr.Face, ip);
176 
177  el.CalcShape(eip, shape);
178 
179  add(elvect, val, shape, elvect);
180  }
181 }
182 
184  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
185 {
186  int dim = el.GetDim()+1;
187  int dof = el.GetDof();
188  Vector nor(dim), Qvec;
189 
190  shape.SetSize(dof);
191  elvect.SetSize(dof);
192  elvect = 0.0;
193 
194  const IntegrationRule *ir = IntRule;
195  if (ir == NULL)
196  {
197  int intorder = oa * el.GetOrder() + ob; // <----------
198  ir = &IntRules.Get(el.GetGeomType(), intorder);
199  }
200 
201  for (int i = 0; i < ir->GetNPoints(); i++)
202  {
203  const IntegrationPoint &ip = ir->IntPoint(i);
204 
205  Tr.SetIntPoint(&ip);
206  if (dim > 1)
207  {
208  CalcOrtho(Tr.Jacobian(), nor);
209  }
210  else
211  {
212  nor[0] = 1.0;
213  }
214  Q.Eval(Qvec, Tr, ip);
215 
216  el.CalcShape(ip, shape);
217 
218  elvect.Add(ip.weight*(Qvec*nor), shape);
219  }
220 }
221 
223  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
224 {
225  int dim = el.GetDim()+1;
226  int dof = el.GetDof();
227  Vector tangent(dim), Qvec;
228 
229  shape.SetSize(dof);
230  elvect.SetSize(dof);
231  elvect = 0.0;
232 
233  if (dim != 2)
234  {
235  mfem_error("These methods make sense only in 2D problems.");
236  }
237 
238  const IntegrationRule *ir = IntRule;
239  if (ir == NULL)
240  {
241  int intorder = oa * el.GetOrder() + ob; // <----------
242  ir = &IntRules.Get(el.GetGeomType(), intorder);
243  }
244 
245  for (int i = 0; i < ir->GetNPoints(); i++)
246  {
247  const IntegrationPoint &ip = ir->IntPoint(i);
248 
249  Tr.SetIntPoint(&ip);
250  const DenseMatrix &Jac = Tr.Jacobian();
251  tangent(0) = Jac(0,0);
252  tangent(1) = Jac(1,0);
253 
254  Q.Eval(Qvec, Tr, ip);
255 
256  el.CalcShape(ip, shape);
257 
258  add(elvect, ip.weight*(Qvec*tangent), shape, elvect);
259  }
260 }
261 
263  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
264 {
265  int vdim = Q.GetVDim();
266  int dof = el.GetDof();
267 
268  double val,cf;
269 
270  shape.SetSize(dof); // vector of size dof
271 
272  elvect.SetSize(dof * vdim);
273  elvect = 0.0;
274 
275  const IntegrationRule *ir = IntRule;
276  if (ir == NULL)
277  {
278  int intorder = 2*el.GetOrder();
279  ir = &IntRules.Get(el.GetGeomType(), intorder);
280  }
281 
282  for (int i = 0; i < ir->GetNPoints(); i++)
283  {
284  const IntegrationPoint &ip = ir->IntPoint(i);
285 
286  Tr.SetIntPoint (&ip);
287  val = Tr.Weight();
288 
289  el.CalcShape(ip, shape);
290  Q.Eval (Qvec, Tr, ip);
291 
292  for (int k = 0; k < vdim; k++)
293  {
294  cf = val * Qvec(k);
295 
296  for (int s = 0; s < dof; s++)
297  {
298  elvect(dof*k+s) += ip.weight * cf * shape(s);
299  }
300  }
301  }
302 }
303 
305  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
306 {
307  MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
308  int vdim = Q.GetVDim();
309  int dof = fe.GetDof();
310 
311  shape.SetSize(dof);
312  fe.CalcPhysShape(Trans, shape);
313 
314  vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
315 
316  elvect.SetSize(dof*vdim);
317  DenseMatrix elvec_as_mat(elvect.GetData(), dof, vdim);
318  MultVWt(shape, Qvec, elvec_as_mat);
319 }
320 
322  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
323 {
324  int vdim = Q.GetVDim();
325  int dof = el.GetDof();
326 
327  shape.SetSize(dof);
328  vec.SetSize(vdim);
329 
330  elvect.SetSize(dof * vdim);
331  elvect = 0.0;
332 
333  const IntegrationRule *ir = IntRule;
334  if (ir == NULL)
335  {
336  int intorder = 2*el.GetOrder();
337  ir = &IntRules.Get(el.GetGeomType(), intorder);
338  }
339 
340  for (int i = 0; i < ir->GetNPoints(); i++)
341  {
342  const IntegrationPoint &ip = ir->IntPoint(i);
343 
344  Q.Eval(vec, Tr, ip);
345  Tr.SetIntPoint (&ip);
346  vec *= Tr.Weight() * ip.weight;
347  el.CalcShape(ip, shape);
348  for (int k = 0; k < vdim; k++)
349  for (int s = 0; s < dof; s++)
350  {
351  elvect(dof*k+s) += vec(k) * shape(s);
352  }
353  }
354 }
355 
357  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
358 {
359  int vdim = Q.GetVDim();
360  int dof = el.GetDof();
361 
362  shape.SetSize(dof);
363  vec.SetSize(vdim);
364 
365  elvect.SetSize(dof * vdim);
366  elvect = 0.0;
367 
368  const IntegrationRule *ir = IntRule;
369  if (ir == NULL)
370  {
371  int intorder = 2*el.GetOrder();
372  ir = &IntRules.Get(Tr.GetGeometryType(), intorder);
373  }
374 
375  for (int i = 0; i < ir->GetNPoints(); i++)
376  {
377  const IntegrationPoint &ip = ir->IntPoint(i);
378 
379  // Set the integration point in the face and the neighboring element
380  Tr.SetAllIntPoints(&ip);
381 
382  // Access the neighboring element's integration point
383  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
384 
385  // Use Tr transformation in case Q depends on boundary attribute
386  Q.Eval(vec, Tr, ip);
387  vec *= Tr.Weight() * ip.weight;
388  el.CalcShape(eip, shape);
389  for (int k = 0; k < vdim; k++)
390  {
391  for (int s = 0; s < dof; s++)
392  {
393  elvect(dof*k+s) += vec(k) * shape(s);
394  }
395  }
396  }
397 }
398 
400  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
401 {
402  int dof = el.GetDof();
403  int spaceDim = Tr.GetSpaceDim();
404  int vdim = std::max(spaceDim, el.GetVDim());
405 
406  vshape.SetSize(dof,vdim);
407  vec.SetSize(vdim);
408 
409  elvect.SetSize(dof);
410  elvect = 0.0;
411 
412  const IntegrationRule *ir = IntRule;
413  if (ir == NULL)
414  {
415  // int intorder = 2*el.GetOrder() - 1; // ok for O(h^{k+1}) conv. in L2
416  int intorder = 2*el.GetOrder();
417  ir = &IntRules.Get(el.GetGeomType(), intorder);
418  }
419 
420  for (int i = 0; i < ir->GetNPoints(); i++)
421  {
422  const IntegrationPoint &ip = ir->IntPoint(i);
423 
424  Tr.SetIntPoint (&ip);
425  el.CalcVShape(Tr, vshape);
426 
427  QF.Eval (vec, Tr, ip);
428  vec *= ip.weight * Tr.Weight();
429  vshape.AddMult (vec, elvect);
430  }
431 }
432 
434  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
435 {
436  MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
437  int dof = fe.GetDof();
438  int spaceDim = Trans.GetSpaceDim();
439 
440  vshape.SetSize(dof, spaceDim);
441  fe.CalcPhysVShape(Trans, vshape);
442 
443  vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
444 
445  elvect.SetSize(dof);
446  vshape.Mult(vec, elvect);
447 }
448 
450  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
451 {
452  int dof = el.GetDof();
453  int spaceDim = Tr.GetSpaceDim();
454  int n=(spaceDim == 3)? spaceDim : 1;
455  curlshape.SetSize(dof,n);
456  vec.SetSize(n);
457 
458  elvect.SetSize(dof);
459  elvect = 0.0;
460 
461  const IntegrationRule *ir = IntRule;
462  if (ir == NULL)
463  {
464  int intorder = 2*el.GetOrder();
465  ir = &IntRules.Get(el.GetGeomType(), intorder);
466  }
467 
468  for (int i = 0; i < ir->GetNPoints(); i++)
469  {
470  const IntegrationPoint &ip = ir->IntPoint(i);
471 
472  Tr.SetIntPoint (&ip);
473  el.CalcPhysCurlShape(Tr, curlshape);
474  QF->Eval(vec, Tr, ip);
475 
476  vec *= ip.weight * Tr.Weight();
477  curlshape.AddMult (vec, elvect);
478  }
479 }
480 
482  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
483 {
484  int spaceDim = Trans.GetSpaceDim();
485  MFEM_ASSERT(vec_delta != NULL,
486  "coefficient must be VectorDeltaCoefficient");
487  int dof = fe.GetDof();
488  int n=(spaceDim == 3)? spaceDim : 1;
489  vec.SetSize(n);
490  curlshape.SetSize(dof, n);
491  elvect.SetSize(dof);
492  fe.CalcPhysCurlShape(Trans, curlshape);
493 
494  vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
495  curlshape.Mult(vec, elvect);
496 }
497 
499  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
500 {
501  int dof = el.GetDof();
502 
503  divshape.SetSize(dof); // vector of size dof
504  elvect.SetSize(dof);
505  elvect = 0.0;
506 
507  const IntegrationRule *ir = IntRule;
508  if (ir == NULL)
509  {
510  int intorder = 2 * el.GetOrder();
511  ir = &IntRules.Get(el.GetGeomType(), intorder);
512  }
513 
514  for (int i = 0; i < ir->GetNPoints(); i++)
515  {
516  const IntegrationPoint &ip = ir->IntPoint(i);
517 
518  Tr.SetIntPoint (&ip);
519  double val = Tr.Weight() * Q.Eval(Tr, ip);
520  el.CalcPhysDivShape(Tr, divshape);
521 
522  add(elvect, ip.weight * val, divshape, elvect);
523  }
524 }
525 
527  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
528 {
529  MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
530  elvect.SetSize(fe.GetDof());
531  fe.CalcPhysDivShape(Trans, elvect);
532  elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
533 }
534 
536  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
537 {
538  int dim = el.GetDim()+1;
539  int dof = el.GetDof();
540 
541  shape.SetSize (dof);
542  nor.SetSize (dim);
543  elvect.SetSize (dim*dof);
544 
545  const IntegrationRule *ir = IntRule;
546  if (ir == NULL)
547  {
548  ir = &IntRules.Get(el.GetGeomType(), el.GetOrder() + 1);
549  }
550 
551  elvect = 0.0;
552  for (int i = 0; i < ir->GetNPoints(); i++)
553  {
554  const IntegrationPoint &ip = ir->IntPoint(i);
555  Tr.SetIntPoint (&ip);
556  CalcOrtho(Tr.Jacobian(), nor);
557  el.CalcShape (ip, shape);
558  nor *= Sign * ip.weight * F -> Eval (Tr, ip);
559  for (int j = 0; j < dof; j++)
560  for (int k = 0; k < dim; k++)
561  {
562  elvect(dof*k+j) += nor(k) * shape(j);
563  }
564  }
565 }
566 
567 
569  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
570 {
571  int dof = el.GetDof();
572 
573  shape.SetSize(dof);
574  elvect.SetSize(dof);
575  elvect = 0.0;
576 
577  const IntegrationRule *ir = IntRule;
578  if (ir == NULL)
579  {
580  int intorder = oa * el.GetOrder() + ob; // <----------
581  ir = &IntRules.Get(el.GetGeomType(), intorder);
582  }
583 
584  for (int i = 0; i < ir->GetNPoints(); i++)
585  {
586  const IntegrationPoint &ip = ir->IntPoint(i);
587  el.CalcShape(ip, shape);
588 
589  double val = ip.weight;
590  if (F)
591  {
592  Tr.SetIntPoint (&ip);
593  val *= F->Eval(Tr, ip);
594  }
595 
596  elvect.Add(val, shape);
597  }
598 }
599 
601  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
602 {
603  int dof = el.GetDof();
604  int dim = el.GetDim();
605  int vdim = el.GetVDim();
606  DenseMatrix vshape(dof, vdim);
607  Vector f_loc(3);
608  Vector f_hat(2);
609 
610  MFEM_VERIFY(vdim == 2, "VectorFEBoundaryTangentLFIntegrator "
611  "must be called with vector basis functions of dimension 2.");
612 
613  elvect.SetSize(dof);
614  elvect = 0.0;
615 
616  const IntegrationRule *ir = IntRule;
617  if (ir == NULL)
618  {
619  int intorder = oa * el.GetOrder() + ob; // <----------
620  ir = &IntRules.Get(el.GetGeomType(), intorder);
621  }
622 
623  for (int i = 0; i < ir->GetNPoints(); i++)
624  {
625  const IntegrationPoint &ip = ir->IntPoint(i);
626 
627  el.CalcVShape(ip, vshape);
628 
629  Tr.SetIntPoint(&ip);
630  f.Eval(f_loc, Tr, ip);
631 
632  if (dim == 2)
633  {
634  Tr.Jacobian().MultTranspose(f_loc, f_hat);
635  }
636  else if (dim == 1)
637  {
638  const DenseMatrix & J = Tr.Jacobian();
639  f_hat(0) = J(0,0) * f_loc(0) + J(1,0) * f_loc(1);
640  f_hat(1) = f_loc(2);
641  }
642  else
643  {
644  f_hat(0) = f_loc(1);
645  f_hat(1) = f_loc(2);
646  }
647 
648  Swap<double>(f_hat(0), f_hat(1));
649  f_hat(0) = -f_hat(0);
650  f_hat *= ip.weight;
651 
652  vshape.AddMult(f_hat, elvect);
653  }
654 }
655 
657  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
658 {
659  mfem_error("BoundaryFlowIntegrator::AssembleRHSElementVect\n"
660  " is not implemented as boundary integrator!\n"
661  " Use LinearForm::AddBdrFaceIntegrator instead of\n"
662  " LinearForm::AddBoundaryIntegrator.");
663 }
664 
666  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
667 {
668  int dim, ndof, order;
669  double un, w, vu_data[3], nor_data[3];
670 
671  dim = el.GetDim();
672  ndof = el.GetDof();
673  Vector vu(vu_data, dim), nor(nor_data, dim);
674 
675  const IntegrationRule *ir = IntRule;
676  if (ir == NULL)
677  {
678  // Assuming order(u)==order(mesh)
679  order = Tr.Elem1->OrderW() + 2*el.GetOrder();
680  if (el.Space() == FunctionSpace::Pk)
681  {
682  order++;
683  }
684  ir = &IntRules.Get(Tr.GetGeometryType(), order);
685  }
686 
687  shape.SetSize(ndof);
688  elvect.SetSize(ndof);
689  elvect = 0.0;
690 
691  for (int p = 0; p < ir->GetNPoints(); p++)
692  {
693  const IntegrationPoint &ip = ir->IntPoint(p);
694 
695  // Set the integration point in the face and the neighboring element
696  Tr.SetAllIntPoints(&ip);
697 
698  // Access the neighboring element's integration point
699  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
700  el.CalcShape(eip, shape);
701 
702  // Use Tr.Elem1 transformation for u so that it matches the coefficient
703  // used with the ConvectionIntegrator and/or the DGTraceIntegrator.
704  u->Eval(vu, *Tr.Elem1, eip);
705 
706  if (dim == 1)
707  {
708  nor(0) = 2*eip.x - 1.0;
709  }
710  else
711  {
712  CalcOrtho(Tr.Jacobian(), nor);
713  }
714 
715  un = vu * nor;
716  w = 0.5*alpha*un - beta*fabs(un);
717  w *= ip.weight*f->Eval(Tr, ip);
718  elvect.Add(w, shape);
719  }
720 }
721 
723  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
724 {
725  mfem_error("DGDirichletLFIntegrator::AssembleRHSElementVect");
726 }
727 
729  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
730 {
731  int dim, ndof;
732  bool kappa_is_nonzero = (kappa != 0.);
733  double w;
734 
735  dim = el.GetDim();
736  ndof = el.GetDof();
737 
738  nor.SetSize(dim);
739  nh.SetSize(dim);
740  ni.SetSize(dim);
741  adjJ.SetSize(dim);
742  if (MQ)
743  {
744  mq.SetSize(dim);
745  }
746 
747  shape.SetSize(ndof);
748  dshape.SetSize(ndof, dim);
749  dshape_dn.SetSize(ndof);
750 
751  elvect.SetSize(ndof);
752  elvect = 0.0;
753 
754  const IntegrationRule *ir = IntRule;
755  if (ir == NULL)
756  {
757  // a simple choice for the integration order; is this OK?
758  int order = 2*el.GetOrder();
759  ir = &IntRules.Get(Tr.GetGeometryType(), order);
760  }
761 
762  for (int p = 0; p < ir->GetNPoints(); p++)
763  {
764  const IntegrationPoint &ip = ir->IntPoint(p);
765 
766  // Set the integration point in the face and the neighboring element
767  Tr.SetAllIntPoints(&ip);
768 
769  // Access the neighboring element's integration point
770  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
771 
772  if (dim == 1)
773  {
774  nor(0) = 2*eip.x - 1.0;
775  }
776  else
777  {
778  CalcOrtho(Tr.Jacobian(), nor);
779  }
780 
781  el.CalcShape(eip, shape);
782  el.CalcDShape(eip, dshape);
783 
784  // compute uD through the face transformation
785  w = ip.weight * uD->Eval(Tr, ip) / Tr.Elem1->Weight();
786  if (!MQ)
787  {
788  if (Q)
789  {
790  w *= Q->Eval(*Tr.Elem1, eip);
791  }
792  ni.Set(w, nor);
793  }
794  else
795  {
796  nh.Set(w, nor);
797  MQ->Eval(mq, *Tr.Elem1, eip);
798  mq.MultTranspose(nh, ni);
799  }
801  adjJ.Mult(ni, nh);
802 
804  elvect.Add(sigma, dshape_dn);
805 
806  if (kappa_is_nonzero)
807  {
808  elvect.Add(kappa*(ni*nor), shape);
809  }
810  }
811 }
812 
814  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
815 {
816  mfem_error("DGElasticityDirichletLFIntegrator::AssembleRHSElementVect");
817 }
818 
820  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
821 {
822  MFEM_ASSERT(Tr.Elem2No < 0, "interior boundary is not supported");
823 
824 #ifdef MFEM_THREAD_SAFE
825  Vector shape;
829  Vector nor;
832  Vector u_dir;
833 #endif
834 
835  const int dim = el.GetDim();
836  const int ndofs = el.GetDof();
837  const int nvdofs = dim*ndofs;
838 
839  elvect.SetSize(nvdofs);
840  elvect = 0.0;
841 
842  adjJ.SetSize(dim);
843  shape.SetSize(ndofs);
844  dshape.SetSize(ndofs, dim);
845  dshape_ps.SetSize(ndofs, dim);
846  nor.SetSize(dim);
847  dshape_dn.SetSize(ndofs);
848  dshape_du.SetSize(ndofs);
849  u_dir.SetSize(dim);
850 
851  const IntegrationRule *ir = IntRule;
852  if (ir == NULL)
853  {
854  const int order = 2*el.GetOrder(); // <-----
855  ir = &IntRules.Get(Tr.GetGeometryType(), order);
856  }
857 
858  for (int pi = 0; pi < ir->GetNPoints(); ++pi)
859  {
860  const IntegrationPoint &ip = ir->IntPoint(pi);
861 
862  // Set the integration point in the face and the neighboring element
863  Tr.SetAllIntPoints(&ip);
864 
865  // Access the neighboring element's integration point
866  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
867 
868  // Evaluate the Dirichlet b.c. using the face transformation.
869  uD.Eval(u_dir, Tr, ip);
870 
871  el.CalcShape(eip, shape);
872  el.CalcDShape(eip, dshape);
873 
875  Mult(dshape, adjJ, dshape_ps);
876 
877  if (dim == 1)
878  {
879  nor(0) = 2*eip.x - 1.0;
880  }
881  else
882  {
883  CalcOrtho(Tr.Jacobian(), nor);
884  }
885 
886  double wL, wM, jcoef;
887  {
888  const double w = ip.weight / Tr.Elem1->Weight();
889  wL = w * lambda->Eval(*Tr.Elem1, eip);
890  wM = w * mu->Eval(*Tr.Elem1, eip);
891  jcoef = kappa * (wL + 2.0*wM) * (nor*nor);
892  dshape_ps.Mult(nor, dshape_dn);
893  dshape_ps.Mult(u_dir, dshape_du);
894  }
895 
896  // alpha < uD, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
897  // + kappa < h^{-1} (lambda + 2 mu) uD, v >
898 
899  // i = idof + ndofs * im
900  // v_phi(i,d) = delta(im,d) phi(idof)
901  // div(v_phi(i)) = dphi(idof,im)
902  // (grad(v_phi(i)))(k,l) = delta(im,k) dphi(idof,l)
903  //
904  // term 1:
905  // alpha < uD, lambda div(v_phi(i)) n >
906  // alpha lambda div(v_phi(i)) (uD.n) =
907  // alpha lambda dphi(idof,im) (uD.n) --> quadrature -->
908  // ip.weight/det(J1) alpha lambda (uD.nor) dshape_ps(idof,im) =
909  // alpha * wL * (u_dir*nor) * dshape_ps(idof,im)
910  // term 2:
911  // < alpha uD, mu grad(v_phi(i)).n > =
912  // alpha mu uD^T grad(v_phi(i)) n =
913  // alpha mu uD(k) delta(im,k) dphi(idof,l) n(l) =
914  // alpha mu uD(im) dphi(idof,l) n(l) --> quadrature -->
915  // ip.weight/det(J1) alpha mu uD(im) dshape_ps(idof,l) nor(l) =
916  // alpha * wM * u_dir(im) * dshape_dn(idof)
917  // term 3:
918  // < alpha uD, mu (grad(v_phi(i)))^T n > =
919  // alpha mu n^T grad(v_phi(i)) uD =
920  // alpha mu n(k) delta(im,k) dphi(idof,l) uD(l) =
921  // alpha mu n(im) dphi(idof,l) uD(l) --> quadrature -->
922  // ip.weight/det(J1) alpha mu nor(im) dshape_ps(idof,l) uD(l) =
923  // alpha * wM * nor(im) * dshape_du(idof)
924  // term j:
925  // < kappa h^{-1} (lambda + 2 mu) uD, v_phi(i) > =
926  // kappa/h (lambda + 2 mu) uD(k) v_phi(i,k) =
927  // kappa/h (lambda + 2 mu) uD(k) delta(im,k) phi(idof) =
928  // kappa/h (lambda + 2 mu) uD(im) phi(idof) --> quadrature -->
929  // [ 1/h = |nor|/det(J1) ]
930  // ip.weight/det(J1) |nor|^2 kappa (lambda + 2 mu) uD(im) phi(idof) =
931  // jcoef * u_dir(im) * shape(idof)
932 
933  wM *= alpha;
934  const double t1 = alpha * wL * (u_dir*nor);
935  for (int im = 0, i = 0; im < dim; ++im)
936  {
937  const double t2 = wM * u_dir(im);
938  const double t3 = wM * nor(im);
939  const double tj = jcoef * u_dir(im);
940  for (int idof = 0; idof < ndofs; ++idof, ++i)
941  {
942  elvect(i) += (t1*dshape_ps(idof,im) + t2*dshape_dn(idof) +
943  t3*dshape_du(idof) + tj*shape(idof));
944  }
945  }
946  }
947 }
948 
950  const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
951 {
952  const IntegrationRule *ir =
954 
955  const int nqp = ir->GetNPoints();
956  const int vdim = vqfc.GetVDim();
957  const int ndofs = fe.GetDof();
958  Vector shape(ndofs);
959  Vector temp(vdim);
960  elvect.SetSize(vdim * ndofs);
961  elvect = 0.0;
962  for (int q = 0; q < nqp; q++)
963  {
964  const IntegrationPoint &ip = ir->IntPoint(q);
965  Tr.SetIntPoint(&ip);
966  const double w = Tr.Weight() * ip.weight;
967  vqfc.Eval(temp, Tr, ip);
968  fe.CalcShape(ip, shape);
969  for (int ind = 0; ind < vdim; ind++)
970  {
971  for (int nd = 0; nd < ndofs; nd++)
972  {
973  elvect(nd + ind * ndofs) += w * shape(nd) * temp(ind);
974  }
975  }
976  }
977 }
978 
981  Vector &elvect)
982 {
983  const IntegrationRule *ir =
985 
986  const int nqp = ir->GetNPoints();
987  const int ndofs = fe.GetDof();
988  Vector shape(ndofs);
989  elvect.SetSize(ndofs);
990  elvect = 0.0;
991  for (int q = 0; q < nqp; q++)
992  {
993  const IntegrationPoint &ip = ir->IntPoint(q);
994  Tr.SetIntPoint (&ip);
995  const double w = Tr.Weight() * ip.weight;
996  double temp = qfc.Eval(Tr, ip);
997  fe.CalcShape(ip, shape);
998  shape *= (w * temp);
999  elvect += shape;
1000  }
1001 }
1002 
1003 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:104
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:568
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:526
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_base.cpp:39
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:2742
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:976
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:498
void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in physical space at the po...
Definition: fe_base.cpp:57
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:656
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:813
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:72
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:63
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2074
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:262
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2288
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:120
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:300
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
Definition: fe_base.cpp:70
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:188
const QuadratureFunction & GetQuadFunction() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
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:154
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition: eltrans.hpp:566
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:321
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...
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:798
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp:55
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:481
MatrixCoefficient * MQ
Definition: lininteg.hpp:448
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Definition: fe_base.cpp:191
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void AddMult(const Vector &x, Vector &y) const
y += A.x
Definition: densemat.cpp:211
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:222
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:722
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
Definition: fe_base.hpp:412
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
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_base.hpp:323
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition: fe_base.cpp:181
DeltaCoefficient * delta
Definition: lininteg.hpp:54
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition: eltrans.hpp:577
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:535
Polynomials of order k.
Definition: fe_base.hpp:223
virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:267
const QuadratureFunction & GetQuadFunction() const
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:600
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:252
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:433
int dim
Definition: ex24.cpp:53
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:949
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:183
ElementTransformation * Elem1
Definition: eltrans.hpp:522
ElementTransformation * Face
Definition: eltrans.hpp:523
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:449
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:979
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:160
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:399
RefCoord s[3]
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:304
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:32
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
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:379