MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
lininteg.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#include "fem.hpp"
13#include <cmath>
14#include "intrules.hpp"
15
16namespace mfem
17{
19 const Array<int> &markers,
20 Vector &b)
21{
22 MFEM_ABORT("Not supported.");
23}
24
26 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
27{
28 mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
29}
30
32 const FiniteElement &el1, const FiniteElement &el2,
34{
35 mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
36}
37
39 : DeltaLFIntegrator(QF), Q(QF), oa(a), ob(b)
40{
41 static Kernels kernels;
42}
43
45 const IntegrationRule *ir)
46 : DeltaLFIntegrator(QF, ir), Q(QF), oa(1), ob(1)
47{
48 static Kernels kernels;
49}
50
53 Vector &elvect)
54{
55 int dof = el.GetDof();
56
57 shape.SetSize(dof); // vector of size dof
58 elvect.SetSize(dof);
59 elvect = 0.0;
60
61 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
62
63 if (ir == NULL)
64 {
65 // ir = &IntRules.Get(el.GetGeomType(),
66 // oa * el.GetOrder() + ob + Tr.OrderW());
67 ir = &IntRules.Get(el.GetGeomType(), oa * el.GetOrder() + ob);
68 }
69
70 for (int i = 0; i < ir->GetNPoints(); i++)
71 {
72 const IntegrationPoint &ip = ir->IntPoint(i);
73
74 Tr.SetIntPoint (&ip);
75 real_t val = Tr.Weight() * Q.Eval(Tr, ip);
76
77 el.CalcPhysShape(Tr, shape);
78
79 add(elvect, ip.weight * val, shape, elvect);
80 }
81}
82
84 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
85{
86 MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
87 elvect.SetSize(fe.GetDof());
88 fe.CalcPhysShape(Trans, elvect);
89 elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
90}
91
93 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
94{
95 int dof = el.GetDof();
96 int spaceDim = Tr.GetSpaceDim();
97
98 dshape.SetSize(dof, spaceDim);
99
100 elvect.SetSize(dof);
101 elvect = 0.0;
102
103 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
104 if (ir == NULL)
105 {
106 int intorder = 2 * el.GetOrder();
107 ir = &IntRules.Get(el.GetGeomType(), intorder);
108 }
109
110 for (int i = 0; i < ir->GetNPoints(); i++)
111 {
112 const IntegrationPoint &ip = ir->IntPoint(i);
113
114 Tr.SetIntPoint(&ip);
115 el.CalcPhysDShape(Tr, dshape);
116
117 Q.Eval(Qvec, Tr, ip);
118 Qvec *= ip.weight * Tr.Weight();
119
120 dshape.AddMult(Qvec, elvect);
121 }
122}
123
125 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
126{
127 MFEM_ASSERT(vec_delta != NULL,"coefficient must be VectorDeltaCoefficient");
128 int dof = fe.GetDof();
129 int spaceDim = Trans.GetSpaceDim();
130
131 dshape.SetSize(dof, spaceDim);
132 fe.CalcPhysDShape(Trans, dshape);
133
134 vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
135
136 elvect.SetSize(dof);
137 dshape.Mult(Qvec, elvect);
138}
139
141 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
142{
143 int dof = el.GetDof();
144
145 shape.SetSize(dof); // vector of size dof
146 elvect.SetSize(dof);
147 elvect = 0.0;
148
149 const IntegrationRule *ir = IntRule;
150 if (ir == NULL)
151 {
152 int intorder = oa * el.GetOrder() + ob; // <----------
153 ir = &IntRules.Get(el.GetGeomType(), intorder);
154 }
155
156 for (int i = 0; i < ir->GetNPoints(); i++)
157 {
158 const IntegrationPoint &ip = ir->IntPoint(i);
159
160 Tr.SetIntPoint (&ip);
161 real_t val = Tr.Weight() * Q.Eval(Tr, ip);
162
163 el.CalcShape(ip, shape);
164
165 add(elvect, ip.weight * val, shape, elvect);
166 }
167}
168
170 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
171{
172 int dof = el.GetDof();
173
174 shape.SetSize(dof); // vector of size dof
175 elvect.SetSize(dof);
176 elvect = 0.0;
177
178 const IntegrationRule *ir = IntRule;
179 if (ir == NULL)
180 {
181 int intorder = oa * el.GetOrder() + ob; // <------ user control
182 ir = &IntRules.Get(Tr.FaceGeom, intorder); // of integration order
183 }
184
185 for (int i = 0; i < ir->GetNPoints(); i++)
186 {
187 const IntegrationPoint &ip = ir->IntPoint(i);
188
189 // Set the integration point in the face and the neighboring element
190 Tr.SetAllIntPoints(&ip);
191
192 // Access the neighboring element's integration point
193 const IntegrationPoint &eip = Tr.GetElement1IntPoint();
194
195 real_t val = Tr.Face->Weight() * ip.weight * Q.Eval(*Tr.Face, ip);
196
197 el.CalcShape(eip, shape);
198
199 add(elvect, val, shape, elvect);
200 }
201}
202
204 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
205{
206 int dim = el.GetDim()+1;
207 int dof = el.GetDof();
208 Vector nor(dim), Qvec;
209
210 shape.SetSize(dof);
211 elvect.SetSize(dof);
212 elvect = 0.0;
213
214 const IntegrationRule *ir = IntRule;
215 if (ir == NULL)
216 {
217 int intorder = oa * el.GetOrder() + ob; // <----------
218 ir = &IntRules.Get(el.GetGeomType(), intorder);
219 }
220
221 for (int i = 0; i < ir->GetNPoints(); i++)
222 {
223 const IntegrationPoint &ip = ir->IntPoint(i);
224
225 Tr.SetIntPoint(&ip);
226 if (dim > 1)
227 {
228 CalcOrtho(Tr.Jacobian(), nor);
229 }
230 else
231 {
232 nor[0] = 1.0;
233 }
234 Q.Eval(Qvec, Tr, ip);
235
236 el.CalcShape(ip, shape);
237
238 elvect.Add(ip.weight*(Qvec*nor), shape);
239 }
240}
241
243 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
244{
245 int dim = el.GetDim()+1;
246 int dof = el.GetDof();
247 Vector tangent(dim), Qvec;
248
249 shape.SetSize(dof);
250 elvect.SetSize(dof);
251 elvect = 0.0;
252
253 if (dim != 2)
254 {
255 mfem_error("These methods make sense only in 2D problems.");
256 }
257
258 const IntegrationRule *ir = IntRule;
259 if (ir == NULL)
260 {
261 int intorder = oa * el.GetOrder() + ob; // <----------
262 ir = &IntRules.Get(el.GetGeomType(), intorder);
263 }
264
265 for (int i = 0; i < ir->GetNPoints(); i++)
266 {
267 const IntegrationPoint &ip = ir->IntPoint(i);
268
269 Tr.SetIntPoint(&ip);
270 const DenseMatrix &Jac = Tr.Jacobian();
271 tangent(0) = Jac(0,0);
272 tangent(1) = Jac(1,0);
273
274 Q.Eval(Qvec, Tr, ip);
275
276 el.CalcShape(ip, shape);
277
278 add(elvect, ip.weight*(Qvec*tangent), shape, elvect);
279 }
280}
281
288
290 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
291{
292 int vdim = Q.GetVDim();
293 int dof = el.GetDof();
294
295 real_t val,cf;
296
297 shape.SetSize(dof); // vector of size dof
298
299 elvect.SetSize(dof * vdim);
300 elvect = 0.0;
301
302 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
303 if (ir == NULL)
304 {
305 int intorder = 2*el.GetOrder();
306 ir = &IntRules.Get(el.GetGeomType(), intorder);
307 }
308
309 for (int i = 0; i < ir->GetNPoints(); i++)
310 {
311 const IntegrationPoint &ip = ir->IntPoint(i);
312
313 Tr.SetIntPoint (&ip);
314 val = Tr.Weight();
315
316 el.CalcPhysShape(Tr, shape);
317 Q.Eval (Qvec, Tr, ip);
318
319 for (int k = 0; k < vdim; k++)
320 {
321 cf = val * Qvec(k);
322
323 for (int s = 0; s < dof; s++)
324 {
325 elvect(dof*k+s) += ip.weight * cf * shape(s);
326 }
327 }
328 }
329}
330
332 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
333{
334 MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
335 int vdim = Q.GetVDim();
336 int dof = fe.GetDof();
337
338 shape.SetSize(dof);
339 fe.CalcPhysShape(Trans, shape);
340
341 vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
342
343 elvect.SetSize(dof*vdim);
344 DenseMatrix elvec_as_mat(elvect.GetData(), dof, vdim);
345 MultVWt(shape, Qvec, elvec_as_mat);
346}
347
349 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
350{
351 const int dim = el.GetDim();
352 const int dof = el.GetDof();
353 const int vdim = Q.GetVDim();
354 const int sdim = Tr.GetSpaceDim();
355
356 dshape.SetSize(dof,sdim);
357
358 elvect.SetSize(dof*(vdim/sdim));
359 elvect = 0.0;
360
361 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
362 if (ir == NULL)
363 {
364 int intorder = 2 * el.GetOrder();
365 ir = &IntRules.Get(el.GetGeomType(), intorder);
366 }
367
368 Vector pelvect(dof);
369 Vector part_x(dim);
370
371 for (int q = 0; q < ir->GetNPoints(); q++)
372 {
373 const IntegrationPoint &ip = ir->IntPoint(q);
374
375 Tr.SetIntPoint(&ip);
376 el.CalcPhysDShape(Tr, dshape);
377
378 Q.Eval(Qvec, Tr, ip);
379 Qvec *= ip.weight * Tr.Weight();
380
381 for (int k = 0; k < vdim/sdim; k++)
382 {
383 for (int d=0; d < sdim; ++d) { part_x(d) = Qvec(k*sdim+d); }
384 dshape.Mult(part_x, pelvect);
385 for (int s = 0; s < dof; ++s) { elvect(s+k*dof) += pelvect(s); }
386 }
387 }
388}
389
392{
393 MFEM_ABORT("Not implemented!");
394}
395
397 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
398{
399 int vdim = Q.GetVDim();
400 int dof = el.GetDof();
401
402 shape.SetSize(dof);
403 vec.SetSize(vdim);
404
405 elvect.SetSize(dof * vdim);
406 elvect = 0.0;
407
408 const IntegrationRule *ir = IntRule;
409 if (ir == NULL)
410 {
411 int intorder = 2*el.GetOrder();
412 ir = &IntRules.Get(el.GetGeomType(), intorder);
413 }
414
415 for (int i = 0; i < ir->GetNPoints(); i++)
416 {
417 const IntegrationPoint &ip = ir->IntPoint(i);
418
419 Tr.SetIntPoint (&ip);
420 Q.Eval(vec, Tr, ip);
421 vec *= Tr.Weight() * ip.weight;
422 el.CalcShape(ip, shape);
423 for (int k = 0; k < vdim; k++)
424 for (int s = 0; s < dof; s++)
425 {
426 elvect(dof*k+s) += vec(k) * shape(s);
427 }
428 }
429}
430
432 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
433{
434 int vdim = Q.GetVDim();
435 int dof = el.GetDof();
436
437 shape.SetSize(dof);
438 vec.SetSize(vdim);
439
440 elvect.SetSize(dof * vdim);
441 elvect = 0.0;
442
443 const IntegrationRule *ir = IntRule;
444 if (ir == NULL)
445 {
446 int intorder = 2*el.GetOrder();
447 ir = &IntRules.Get(Tr.GetGeometryType(), intorder);
448 }
449
450 for (int i = 0; i < ir->GetNPoints(); i++)
451 {
452 const IntegrationPoint &ip = ir->IntPoint(i);
453
454 // Set the integration point in the face and the neighboring element
455 Tr.SetAllIntPoints(&ip);
456
457 // Access the neighboring element's integration point
458 const IntegrationPoint &eip = Tr.GetElement1IntPoint();
459
460 // Use Tr transformation in case Q depends on boundary attribute
461 Q.Eval(vec, Tr, ip);
462 vec *= Tr.Weight() * ip.weight;
463 el.CalcShape(eip, shape);
464 for (int k = 0; k < vdim; k++)
465 {
466 for (int s = 0; s < dof; s++)
467 {
468 elvect(dof*k+s) += vec(k) * shape(s);
469 }
470 }
471 }
472}
473
475 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
476{
477 int dof = el.GetDof();
478 int spaceDim = Tr.GetSpaceDim();
479 int vdim = std::max(spaceDim, el.GetRangeDim());
480
481 vshape.SetSize(dof,vdim);
482 vec.SetSize(vdim);
483
484 elvect.SetSize(dof);
485 elvect = 0.0;
486
487 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
488 if (ir == NULL)
489 {
490 // int intorder = 2*el.GetOrder() - 1; // ok for O(h^{k+1}) conv. in L2
491 int intorder = 2*el.GetOrder();
492 ir = &IntRules.Get(el.GetGeomType(), intorder);
493 }
494
495 for (int i = 0; i < ir->GetNPoints(); i++)
496 {
497 const IntegrationPoint &ip = ir->IntPoint(i);
498
499 Tr.SetIntPoint (&ip);
500 el.CalcVShape(Tr, vshape);
501
502 QF.Eval (vec, Tr, ip);
503 vec *= ip.weight * Tr.Weight();
504 vshape.AddMult (vec, elvect);
505 }
506}
507
509 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
510{
511 MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
512 int dof = fe.GetDof();
513 int spaceDim = Trans.GetSpaceDim();
514
515 vshape.SetSize(dof, spaceDim);
516 fe.CalcPhysVShape(Trans, vshape);
517
518 vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
519
520 elvect.SetSize(dof);
521 vshape.Mult(vec, elvect);
522}
523
525 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
526{
527 int dof = el.GetDof();
528 int spaceDim = Tr.GetSpaceDim();
529 int n=(spaceDim == 3)? spaceDim : 1;
530 curlshape.SetSize(dof,n);
531 vec.SetSize(n);
532
533 elvect.SetSize(dof);
534 elvect = 0.0;
535
536 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
537 if (ir == NULL)
538 {
539 int intorder = 2*el.GetOrder();
540 ir = &IntRules.Get(el.GetGeomType(), intorder);
541 }
542
543 for (int i = 0; i < ir->GetNPoints(); i++)
544 {
545 const IntegrationPoint &ip = ir->IntPoint(i);
546
547 Tr.SetIntPoint (&ip);
548 el.CalcPhysCurlShape(Tr, curlshape);
549 QF->Eval(vec, Tr, ip);
550
551 vec *= ip.weight * Tr.Weight();
552 curlshape.AddMult (vec, elvect);
553 }
554}
555
557 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
558{
559 int spaceDim = Trans.GetSpaceDim();
560 MFEM_ASSERT(vec_delta != NULL,
561 "coefficient must be VectorDeltaCoefficient");
562 int dof = fe.GetDof();
563 int n=(spaceDim == 3)? spaceDim : 1;
564 vec.SetSize(n);
565 curlshape.SetSize(dof, n);
566 elvect.SetSize(dof);
567 fe.CalcPhysCurlShape(Trans, curlshape);
568
569 vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
570 curlshape.Mult(vec, elvect);
571}
572
574 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
575{
576 int dof = el.GetDof();
577
578 divshape.SetSize(dof); // vector of size dof
579 elvect.SetSize(dof);
580 elvect = 0.0;
581
582 const IntegrationRule *ir = GetIntegrationRule(el, Tr);
583 if (ir == NULL)
584 {
585 int intorder = 2 * el.GetOrder();
586 ir = &IntRules.Get(el.GetGeomType(), intorder);
587 }
588
589 for (int i = 0; i < ir->GetNPoints(); i++)
590 {
591 const IntegrationPoint &ip = ir->IntPoint(i);
592
593 Tr.SetIntPoint (&ip);
594 real_t val = Tr.Weight() * Q.Eval(Tr, ip);
595 el.CalcPhysDivShape(Tr, divshape);
596
597 add(elvect, ip.weight * val, divshape, elvect);
598 }
599}
600
602 const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
603{
604 MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
605 elvect.SetSize(fe.GetDof());
606 fe.CalcPhysDivShape(Trans, elvect);
607 elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
608}
609
611 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
612{
613 int dim = el.GetDim()+1;
614 int dof = el.GetDof();
615
616 shape.SetSize (dof);
617 nor.SetSize (dim);
618 elvect.SetSize (dim*dof);
619
620 const IntegrationRule *ir = IntRule;
621 if (ir == NULL)
622 {
623 ir = &IntRules.Get(el.GetGeomType(), el.GetOrder() + 1);
624 }
625
626 elvect = 0.0;
627 for (int i = 0; i < ir->GetNPoints(); i++)
628 {
629 const IntegrationPoint &ip = ir->IntPoint(i);
630 Tr.SetIntPoint (&ip);
631 CalcOrtho(Tr.Jacobian(), nor);
632 el.CalcShape (ip, shape);
633 nor *= Sign * ip.weight * F -> Eval (Tr, ip);
634 for (int j = 0; j < dof; j++)
635 for (int k = 0; k < dim; k++)
636 {
637 elvect(dof*k+j) += nor(k) * shape(j);
638 }
639 }
640}
641
642
644 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
645{
646 int dof = el.GetDof();
647
648 shape.SetSize(dof);
649 elvect.SetSize(dof);
650 elvect = 0.0;
651
652 const IntegrationRule *ir = IntRule;
653 if (ir == NULL)
654 {
655 int intorder = oa * el.GetOrder() + ob; // <----------
656 ir = &IntRules.Get(el.GetGeomType(), intorder);
657 }
658
659 for (int i = 0; i < ir->GetNPoints(); i++)
660 {
661 const IntegrationPoint &ip = ir->IntPoint(i);
662 el.CalcShape(ip, shape);
663
664 real_t val = ip.weight;
665 if (F)
666 {
667 Tr.SetIntPoint (&ip);
668 val *= F->Eval(Tr, ip);
669 }
670
671 elvect.Add(val, shape);
672 }
673}
674
676 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
677{
678 int dim = el.GetDim()+1;
679 int dof = el.GetDof();
680 Vector nor(dim), Fvec(dim);
681
682 shape.SetSize(dof);
683 elvect.SetSize(dof);
684 elvect = 0.0;
685
686 const IntegrationRule *ir = IntRule;
687 if (ir == NULL)
688 {
689 int intorder = 2 * el.GetOrder() + Tr.OrderW(); // <----------
690 ir = &IntRules.Get(el.GetGeomType(), intorder);
691 }
692
693 for (int i = 0; i < ir->GetNPoints(); i++)
694 {
695 const IntegrationPoint &ip = ir->IntPoint(i);
696
697 Tr.SetIntPoint(&ip);
698 CalcOrtho(Tr.Jacobian(), nor);
699 F.Eval(Fvec, Tr, ip);
700 real_t val = ip.weight * (Fvec*nor) / Tr.Weight();
701
702 el.CalcShape(ip, shape);
703
704 elvect.Add(val, shape);
705 }
706}
707
709 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
710{
711 int dof = el.GetDof();
712 int dim = el.GetDim();
713 int vdim = el.GetRangeDim();
714 DenseMatrix vshape(dof, vdim);
715 Vector f_loc(3);
716 Vector f_hat(2);
717
718 MFEM_VERIFY(vdim == 2, "VectorFEBoundaryTangentLFIntegrator "
719 "must be called with vector basis functions of dimension 2.");
720
721 elvect.SetSize(dof);
722 elvect = 0.0;
723
724 const IntegrationRule *ir = IntRule;
725 if (ir == NULL)
726 {
727 int intorder = oa * el.GetOrder() + ob; // <----------
728 ir = &IntRules.Get(el.GetGeomType(), intorder);
729 }
730
731 for (int i = 0; i < ir->GetNPoints(); i++)
732 {
733 const IntegrationPoint &ip = ir->IntPoint(i);
734
735 el.CalcVShape(ip, vshape);
736
737 Tr.SetIntPoint(&ip);
738 f.Eval(f_loc, Tr, ip);
739
740 if (dim == 2)
741 {
742 Tr.Jacobian().MultTranspose(f_loc, f_hat);
743 }
744 else if (dim == 1)
745 {
746 const DenseMatrix & J = Tr.Jacobian();
747 f_hat(0) = J(0,0) * f_loc(0) + J(1,0) * f_loc(1);
748 f_hat(1) = f_loc(2);
749 }
750 else
751 {
752 f_hat(0) = f_loc(1);
753 f_hat(1) = f_loc(2);
754 }
755
756 Swap<real_t>(f_hat(0), f_hat(1));
757 f_hat(0) = -f_hat(0);
758 f_hat *= ip.weight;
759
760 vshape.AddMult(f_hat, elvect);
761 }
762}
763
765 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
766{
767 mfem_error("BoundaryFlowIntegrator::AssembleRHSElementVect\n"
768 " is not implemented as boundary integrator!\n"
769 " Use LinearForm::AddBdrFaceIntegrator instead of\n"
770 " LinearForm::AddBoundaryIntegrator.");
771}
772
774 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
775{
776 int dim, ndof, order;
777 real_t un, w, vu_data[3], nor_data[3];
778
779 dim = el.GetDim();
780 ndof = el.GetDof();
781 Vector vu(vu_data, dim), nor(nor_data, dim);
782
783 const IntegrationRule *ir = IntRule;
784 if (ir == NULL)
785 {
786 // Assuming order(u)==order(mesh)
787 order = Tr.Elem1->OrderW() + 2*el.GetOrder();
788 if (el.Space() == FunctionSpace::Pk)
789 {
790 order++;
791 }
792 ir = &IntRules.Get(Tr.GetGeometryType(), order);
793 }
794
795 shape.SetSize(ndof);
796 elvect.SetSize(ndof);
797 elvect = 0.0;
798
799 for (int p = 0; p < ir->GetNPoints(); p++)
800 {
801 const IntegrationPoint &ip = ir->IntPoint(p);
802
803 // Set the integration point in the face and the neighboring element
804 Tr.SetAllIntPoints(&ip);
805
806 // Access the neighboring element's integration point
807 const IntegrationPoint &eip = Tr.GetElement1IntPoint();
808 el.CalcShape(eip, shape);
809
810 // Use Tr.Elem1 transformation for u so that it matches the coefficient
811 // used with the ConvectionIntegrator and/or the DGTraceIntegrator.
812 u->Eval(vu, *Tr.Elem1, eip);
813
814 if (dim == 1)
815 {
816 nor(0) = 2*eip.x - 1.0;
817 }
818 else
819 {
820 CalcOrtho(Tr.Jacobian(), nor);
821 }
822
823 un = vu * nor;
824 w = 0.5*alpha*un - beta*fabs(un);
825 w *= ip.weight*f->Eval(Tr, ip);
826 elvect.Add(w, shape);
827 }
828}
829
831 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
832{
833 mfem_error("DGDirichletLFIntegrator::AssembleRHSElementVect");
834}
835
837 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
838{
839 int dim, ndof;
840 bool kappa_is_nonzero = (kappa != 0.);
841 real_t w;
842
843 dim = el.GetDim();
844 ndof = el.GetDof();
845
846 nor.SetSize(dim);
847 nh.SetSize(dim);
848 ni.SetSize(dim);
850 if (MQ)
851 {
852 mq.SetSize(dim);
853 }
854
855 shape.SetSize(ndof);
856 dshape.SetSize(ndof, dim);
857 dshape_dn.SetSize(ndof);
858
859 elvect.SetSize(ndof);
860 elvect = 0.0;
861
862 const IntegrationRule *ir = IntRule;
863 if (ir == NULL)
864 {
865 // a simple choice for the integration order; is this OK?
866 int order = 2*el.GetOrder();
867 ir = &IntRules.Get(Tr.GetGeometryType(), order);
868 }
869
870 for (int p = 0; p < ir->GetNPoints(); p++)
871 {
872 const IntegrationPoint &ip = ir->IntPoint(p);
873
874 // Set the integration point in the face and the neighboring element
875 Tr.SetAllIntPoints(&ip);
876
877 // Access the neighboring element's integration point
878 const IntegrationPoint &eip = Tr.GetElement1IntPoint();
879
880 if (dim == 1)
881 {
882 nor(0) = 2*eip.x - 1.0;
883 }
884 else
885 {
886 CalcOrtho(Tr.Jacobian(), nor);
887 }
888
889 el.CalcShape(eip, shape);
890 el.CalcDShape(eip, dshape);
891
892 // compute uD through the face transformation
893 w = ip.weight * uD->Eval(Tr, ip) / Tr.Elem1->Weight();
894 if (!MQ)
895 {
896 if (Q)
897 {
898 w *= Q->Eval(*Tr.Elem1, eip);
899 }
900 ni.Set(w, nor);
901 }
902 else
903 {
904 nh.Set(w, nor);
905 MQ->Eval(mq, *Tr.Elem1, eip);
907 }
909 adjJ.Mult(ni, nh);
910
912 elvect.Add(sigma, dshape_dn);
913
914 if (kappa_is_nonzero)
915 {
916 elvect.Add(kappa*(ni*nor), shape);
917 }
918 }
919}
920
922 const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
923{
924 mfem_error("DGElasticityDirichletLFIntegrator::AssembleRHSElementVect");
925}
926
928 const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
929{
930 MFEM_ASSERT(Tr.Elem2No < 0, "interior boundary is not supported");
931
932#ifdef MFEM_THREAD_SAFE
937 Vector nor;
941#endif
942
943 const int dim = el.GetDim();
944 const int ndofs = el.GetDof();
945 const int nvdofs = dim*ndofs;
946
947 elvect.SetSize(nvdofs);
948 elvect = 0.0;
949
951 shape.SetSize(ndofs);
952 dshape.SetSize(ndofs, dim);
953 dshape_ps.SetSize(ndofs, dim);
954 nor.SetSize(dim);
955 dshape_dn.SetSize(ndofs);
956 dshape_du.SetSize(ndofs);
958
959 const IntegrationRule *ir = IntRule;
960 if (ir == NULL)
961 {
962 const int order = 2*el.GetOrder(); // <-----
963 ir = &IntRules.Get(Tr.GetGeometryType(), order);
964 }
965
966 for (int pi = 0; pi < ir->GetNPoints(); ++pi)
967 {
968 const IntegrationPoint &ip = ir->IntPoint(pi);
969
970 // Set the integration point in the face and the neighboring element
971 Tr.SetAllIntPoints(&ip);
972
973 // Access the neighboring element's integration point
974 const IntegrationPoint &eip = Tr.GetElement1IntPoint();
975
976 // Evaluate the Dirichlet b.c. using the face transformation.
977 uD.Eval(u_dir, Tr, ip);
978
979 el.CalcShape(eip, shape);
980 el.CalcDShape(eip, dshape);
981
984
985 if (dim == 1)
986 {
987 nor(0) = 2*eip.x - 1.0;
988 }
989 else
990 {
991 CalcOrtho(Tr.Jacobian(), nor);
992 }
993
994 real_t wL, wM, jcoef;
995 {
996 const real_t w = ip.weight / Tr.Elem1->Weight();
997 wL = w * lambda->Eval(*Tr.Elem1, eip);
998 wM = w * mu->Eval(*Tr.Elem1, eip);
999 jcoef = kappa * (wL + 2.0*wM) * (nor*nor);
1002 }
1003
1004 // alpha < uD, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
1005 // + kappa < h^{-1} (lambda + 2 mu) uD, v >
1006
1007 // i = idof + ndofs * im
1008 // v_phi(i,d) = delta(im,d) phi(idof)
1009 // div(v_phi(i)) = dphi(idof,im)
1010 // (grad(v_phi(i)))(k,l) = delta(im,k) dphi(idof,l)
1011 //
1012 // term 1:
1013 // alpha < uD, lambda div(v_phi(i)) n >
1014 // alpha lambda div(v_phi(i)) (uD.n) =
1015 // alpha lambda dphi(idof,im) (uD.n) --> quadrature -->
1016 // ip.weight/det(J1) alpha lambda (uD.nor) dshape_ps(idof,im) =
1017 // alpha * wL * (u_dir*nor) * dshape_ps(idof,im)
1018 // term 2:
1019 // < alpha uD, mu grad(v_phi(i)).n > =
1020 // alpha mu uD^T grad(v_phi(i)) n =
1021 // alpha mu uD(k) delta(im,k) dphi(idof,l) n(l) =
1022 // alpha mu uD(im) dphi(idof,l) n(l) --> quadrature -->
1023 // ip.weight/det(J1) alpha mu uD(im) dshape_ps(idof,l) nor(l) =
1024 // alpha * wM * u_dir(im) * dshape_dn(idof)
1025 // term 3:
1026 // < alpha uD, mu (grad(v_phi(i)))^T n > =
1027 // alpha mu n^T grad(v_phi(i)) uD =
1028 // alpha mu n(k) delta(im,k) dphi(idof,l) uD(l) =
1029 // alpha mu n(im) dphi(idof,l) uD(l) --> quadrature -->
1030 // ip.weight/det(J1) alpha mu nor(im) dshape_ps(idof,l) uD(l) =
1031 // alpha * wM * nor(im) * dshape_du(idof)
1032 // term j:
1033 // < kappa h^{-1} (lambda + 2 mu) uD, v_phi(i) > =
1034 // kappa/h (lambda + 2 mu) uD(k) v_phi(i,k) =
1035 // kappa/h (lambda + 2 mu) uD(k) delta(im,k) phi(idof) =
1036 // kappa/h (lambda + 2 mu) uD(im) phi(idof) --> quadrature -->
1037 // [ 1/h = |nor|/det(J1) ]
1038 // ip.weight/det(J1) |nor|^2 kappa (lambda + 2 mu) uD(im) phi(idof) =
1039 // jcoef * u_dir(im) * shape(idof)
1040
1041 wM *= alpha;
1042 const real_t t1 = alpha * wL * (u_dir*nor);
1043 for (int im = 0, i = 0; im < dim; ++im)
1044 {
1045 const real_t t2 = wM * u_dir(im);
1046 const real_t t3 = wM * nor(im);
1047 const real_t tj = jcoef * u_dir(im);
1048 for (int idof = 0; idof < ndofs; ++idof, ++i)
1049 {
1050 elvect(i) += (t1*dshape_ps(idof,im) + t2*dshape_dn(idof) +
1051 t3*dshape_du(idof) + tj*shape(idof));
1052 }
1053 }
1054 }
1055}
1056
1057
1058
1060(const FiniteElement &el,
1062 Vector &elvect)
1063{
1064 int n = el.GetDof();
1065 elvect.SetSize(n);
1066 for (int i = 0; i < n; i++)
1067 {
1068 elvect(i) = dist(generator);
1069 }
1070
1071 int iel = Tr.ElementNo;
1072
1073 if (!save_factors || !L[iel])
1074 {
1075 DenseMatrix *M, m;
1076 if (save_factors)
1077 {
1078 L[iel]=new DenseMatrix;
1079 M = L[iel];
1080 }
1081 else
1082 {
1083 M = &m;
1084 }
1085 massinteg.AssembleElementMatrix(el, Tr, *M);
1086 CholeskyFactors chol(M->Data());
1087 chol.Factor(M->Height());
1088 chol.LMult(n,1,elvect.GetData());
1089 }
1090 else
1091 {
1092 CholeskyFactors chol(L[iel]->Data());
1093 chol.LMult(n,1,elvect.GetData());
1094 }
1095}
1096
1097
1099 const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
1100{
1101 const IntegrationRule *ir =
1103
1104 const int nqp = ir->GetNPoints();
1105 const int vdim = vqfc.GetVDim();
1106 const int ndofs = fe.GetDof();
1107 Vector shape(ndofs);
1108 Vector temp(vdim);
1109 elvect.SetSize(vdim * ndofs);
1110 elvect = 0.0;
1111 for (int q = 0; q < nqp; q++)
1112 {
1113 const IntegrationPoint &ip = ir->IntPoint(q);
1114 Tr.SetIntPoint(&ip);
1115 const real_t w = Tr.Weight() * ip.weight;
1116 vqfc.Eval(temp, Tr, ip);
1117 fe.CalcShape(ip, shape);
1118 for (int ind = 0; ind < vdim; ind++)
1119 {
1120 for (int nd = 0; nd < ndofs; nd++)
1121 {
1122 elvect(nd + ind * ndofs) += w * shape(nd) * temp(ind);
1123 }
1124 }
1125 }
1126}
1127
1128
1131 Vector &elvect)
1132{
1133 const IntegrationRule *ir =
1135
1136 const int nqp = ir->GetNPoints();
1137 const int ndofs = fe.GetDof();
1138 Vector shape(ndofs);
1139 elvect.SetSize(ndofs);
1140 elvect = 0.0;
1141 for (int q = 0; q < nqp; q++)
1142 {
1143 const IntegrationPoint &ip = ir->IntPoint(q);
1144 Tr.SetIntPoint (&ip);
1145 const real_t w = Tr.Weight() * ip.weight;
1146 real_t temp = qfc.Eval(Tr, ip);
1147 fe.CalcShape(ip, shape);
1148 shape *= (w * temp);
1149 elvect += shape;
1150 }
1151}
1152
1153}
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:764
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:140
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:203
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:242
bool Factor(int m, real_t TOL=0.0) override
Compute the Cholesky factorization of the current matrix.
void LMult(int m, int n, real_t *X) const
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:830
MatrixCoefficient * MQ
Definition lininteg.hpp:577
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:921
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
Abstract class for integrators that support delta coefficients.
Definition lininteg.hpp:62
VectorDeltaCoefficient * vec_delta
Definition lininteg.hpp:65
DeltaCoefficient * delta
Definition lininteg.hpp:64
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:108
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:158
real_t * Data() const
Returns the matrix data array. Warning: this method casts away constness.
Definition densemat.hpp:122
void AddMult(const Vector &x, Vector &y, const real_t a=1.0) const override
y += a * A.x
Definition densemat.cpp:194
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:116
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:124
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:92
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:51
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Definition lininteg.cpp:38
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:83
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:175
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition eltrans.hpp:111
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:144
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint,...
Definition eltrans.hpp:132
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:750
ElementTransformation * Elem1
Definition eltrans.hpp:791
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition eltrans.hpp:846
ElementTransformation * Face
Definition eltrans.hpp:792
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition eltrans.hpp:835
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
Abstract class for all finite elements.
Definition fe_base.hpp:247
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:50
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
Definition fe_base.hpp:328
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:341
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:202
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:324
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...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:334
int Space() const
Returns the type of FunctionSpace on the element.
Definition fe_base.hpp:351
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:68
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
Definition fe_base.hpp:450
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...
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:81
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:337
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:192
@ Pk
Polynomials of order k.
Definition fe_base.hpp:233
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
const IntegrationRule * GetIntegrationRule() const
Equivalent to GetIntRule, but retained for backward compatibility with applications.
const IntegrationRule * IntRule
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Definition lininteg.cpp:18
void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) override
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 ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient in the element described by T at the point ip.
const QuadratureFunction & GetQuadFunction() const
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition qfunction.hpp:94
void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect) override
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition qspace.hpp:125
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:610
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:396
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
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 ...
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:348
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:390
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:331
VectorDomainLFIntegrator(VectorCoefficient &QF, const IntegrationRule *ir=nullptr)
Constructs a domain integrator with a given VectorCoefficient.
Definition lininteg.cpp:282
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:289
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:643
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:675
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:708
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:556
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:524
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:601
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:573
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition lininteg.cpp:474
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition lininteg.cpp:508
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
const QuadratureFunction & GetQuadFunction() const
void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect) override
Vector data type.
Definition vector.hpp:82
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:341
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:243
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:326
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void CalcOrtho(const DenseMatrix &J, Vector &n)
void mfem_error(const char *msg)
Definition error.cpp:154
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:505
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:414
void Swap(T &a, T &b)
Swap objects of type T. The operation is performed using the most specialized swap function from the ...
Definition array.hpp:738
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
float real_t
Definition config.hpp:46
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
real_t p(const Vector &x, real_t t)