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