MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
coefficient.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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// Implementation of Coefficient class
13
14#include "fem.hpp"
15
16#include <cmath>
17#include <limits>
18
19namespace mfem
20{
21
22using namespace std;
23
24// Given an ElementTransformation and IntegrationPoint in a refined mesh,
25// return the ElementTransformation of the parent coarse element, and set
26// coarse_ip to the location of the original ip within the coarse element.
28 Mesh &coarse_mesh, const ElementTransformation &T,
29 const IntegrationPoint &ip, IntegrationPoint &coarse_ip)
30{
31 const Mesh &fine_mesh = *T.mesh;
32 // Get the element transformation of the coarse element containing the
33 // fine element.
34 int fine_element = T.ElementNo;
36 int coarse_element = cf.embeddings[fine_element].parent;
37 ElementTransformation *coarse_T = coarse_mesh.GetElementTransformation(
38 coarse_element);
39 // Transform the integration point from fine element coordinates to coarse
40 // element coordinates.
42 IntegrationPointTransformation fine_to_coarse;
43 IsoparametricTransformation &emb_tr = fine_to_coarse.Transf;
44 emb_tr.SetIdentityTransformation(geom);
45 emb_tr.SetPointMat(cf.point_matrices[geom](cf.embeddings[fine_element].matrix));
46 fine_to_coarse.Transform(ip, coarse_ip);
47 coarse_T->SetIntPoint(&coarse_ip);
48 return coarse_T;
49}
50
52{
53 QuadratureSpaceBase &qspace = *qf.GetSpace();
54 const int ne = qspace.GetNE();
55 Vector values;
56 for (int iel = 0; iel < ne; ++iel)
57 {
58 qf.GetValues(iel, values);
59 const IntegrationRule &ir = qspace.GetIntRule(iel);
61 for (int iq = 0; iq < ir.Size(); ++iq)
62 {
63 const IntegrationPoint &ip = ir[iq];
64 T.SetIntPoint(&ip);
65 const int iq_p = qspace.GetPermutedIndex(iel, iq);
66 values[iq_p] = Eval(T, ip);
67 }
68 }
69}
70
75
77 const IntegrationPoint & ip)
78{
79 int att = T.Attribute;
80 return (constants(att-1));
81}
82
83void PWCoefficient::InitMap(const Array<int> & attr,
84 const Array<Coefficient*> & coefs)
85{
86 MFEM_VERIFY(attr.Size() == coefs.Size(),
87 "PWCoefficient: "
88 "Attribute and coefficient arrays have incompatible "
89 "dimensions.");
90
91 for (int i=0; i<attr.Size(); i++)
92 {
93 if (coefs[i] != NULL)
94 {
95 UpdateCoefficient(attr[i], *coefs[i]);
96 }
97 }
98}
99
101{
103
104 std::map<int, Coefficient*>::iterator p = pieces.begin();
105 for (; p != pieces.end(); p++)
106 {
107 if (p->second != NULL)
108 {
109 p->second->SetTime(t);
110 }
111 }
112}
113
115 const IntegrationPoint &ip)
116{
117 const int att = T.Attribute;
118 std::map<int, Coefficient*>::const_iterator p = pieces.find(att);
119 if (p != pieces.end())
120 {
121 if ( p->second != NULL)
122 {
123 return p->second->Eval(T, ip);
124 }
125 }
126 return 0.0;
127}
128
130 const IntegrationPoint & ip)
131{
132 real_t x[3];
133 Vector transip(x, 3);
134
135 T.Transform(ip, transip);
136
137 if (Function)
138 {
139 return Function(transip);
140 }
141 else
142 {
143 return TDFunction(transip, GetTime());
144 }
145}
146
153
155 const IntegrationPoint & ip)
156{
157 T.Transform(ip, transip);
158 return sqrt(transip[0] * transip[0] + transip[1] * transip[1]);
159}
160
162 const IntegrationPoint & ip)
163{
164 T.Transform(ip, transip);
165 return atan2(transip[1], transip[0]);
166}
167
169 const IntegrationPoint & ip)
170{
171 T.Transform(ip, transip);
172 return sqrt(transip * transip);
173}
174
176 const IntegrationPoint & ip)
177{
178 T.Transform(ip, transip);
179 return atan2(transip[1], transip[0]);
180}
181
183 const IntegrationPoint & ip)
184{
185 T.Transform(ip, transip);
186 return atan2(sqrt(transip[0] * transip[0] + transip[1] * transip[1]),
187 transip[2]);
188}
189
191 const IntegrationPoint &ip)
192{
193 Mesh *gf_mesh = GridF->FESpace()->GetMesh();
194 if (T.mesh->GetNE() == gf_mesh->GetNE())
195 {
196 return GridF->GetValue(T, ip, Component);
197 }
198 else
199 {
200 IntegrationPoint coarse_ip;
201 ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
202 return GridF->GetValue(*coarse_T, coarse_ip, Component);
203 }
204}
205
210
212{
213 if (Q1) { Q1->SetTime(t); }
214 if (Q2) { Q2->SetTime(t); }
215 this->Coefficient::SetTime(t);
216}
217
219 const IntegrationPoint &ip)
220{
221 if (Q2)
222 {
223 return Transform2(Q1->Eval(T, ip, GetTime()),
224 Q2->Eval(T, ip, GetTime()));
225 }
226 else
227 {
228 return Transform1(Q1->Eval(T, ip, GetTime()));
229 }
230}
231
233{
234 if (weight) { weight->SetTime(t); }
235 this->Coefficient::SetTime(t);
236}
237
239{
240 MFEM_VERIFY(vcenter.Size() <= 3,
241 "SetDeltaCenter::Maximum number of dim supported is 3")
242 for (int i = 0; i < vcenter.Size(); i++) { center[i] = vcenter[i]; }
243 sdim = vcenter.Size();
244}
245
247{
248 vcenter.SetSize(sdim);
249 vcenter = center;
250}
251
253 const IntegrationPoint &ip)
254{
255 real_t w = Scale();
256 return weight ? weight->Eval(T, ip, GetTime())*w : w;
257}
258
260{
261 if (c) { c->SetTime(t); }
262 this->Coefficient::SetTime(t);
263}
264
266 const IntegrationRule &ir)
267{
268 Vector Mi;
269 M.SetSize(vdim, ir.GetNPoints());
270 for (int i = 0; i < ir.GetNPoints(); i++)
271 {
272 M.GetColumnReference(i, Mi);
273 const IntegrationPoint &ip = ir.IntPoint(i);
274 T.SetIntPoint(&ip);
275 Eval(Mi, T, ip);
276 }
277}
278
280{
281 MFEM_VERIFY(vdim == qf.GetVDim(), "Wrong sizes.");
282 QuadratureSpaceBase &qspace = *qf.GetSpace();
283 const int ne = qspace.GetNE();
284 DenseMatrix values;
285 Vector col;
286 for (int iel = 0; iel < ne; ++iel)
287 {
288 qf.GetValues(iel, values);
289 const IntegrationRule &ir = qspace.GetIntRule(iel);
291 for (int iq = 0; iq < ir.Size(); ++iq)
292 {
293 const IntegrationPoint &ip = ir[iq];
294 T.SetIntPoint(&ip);
295 const int iq_p = qspace.GetPermutedIndex(iel, iq);
296 values.GetColumnReference(iq_p, col);
297 Eval(col, T, ip);
298 }
299 }
300}
301
302void PWVectorCoefficient::InitMap(const Array<int> & attr,
303 const Array<VectorCoefficient*> & coefs)
304{
305 MFEM_VERIFY(attr.Size() == coefs.Size(),
306 "PWVectorCoefficient: "
307 "Attribute and coefficient arrays have incompatible "
308 "dimensions.");
309
310 for (int i=0; i<attr.Size(); i++)
311 {
312 if (coefs[i] != NULL)
313 {
314 UpdateCoefficient(attr[i], *coefs[i]);
315 }
316 }
317}
318
320{
321 MFEM_VERIFY(coef.GetVDim() == vdim,
322 "PWVectorCoefficient::UpdateCoefficient: "
323 "VectorCoefficient has incompatible dimension.");
324 pieces[attr] = &coef;
325}
326
328{
330
331 std::map<int, VectorCoefficient*>::iterator p = pieces.begin();
332 for (; p != pieces.end(); p++)
333 {
334 if (p->second != NULL)
335 {
336 p->second->SetTime(t);
337 }
338 }
339}
340
342 const IntegrationPoint &ip)
343{
344 const int att = T.Attribute;
345 std::map<int, VectorCoefficient*>::const_iterator p = pieces.find(att);
346 if (p != pieces.end())
347 {
348 if ( p->second != NULL)
349 {
350 p->second->Eval(V, T, ip);
351 return;
352 }
353 }
354
355 V.SetSize(vdim);
356 V = 0.0;
357}
358
360 const IntegrationPoint &ip)
361{
362 V.SetSize(vdim);
363 T.Transform(ip, V);
364}
365
367 const IntegrationPoint &ip)
368{
369 real_t x[3];
370 Vector transip(x, 3);
371
372 T.Transform(ip, transip);
373
374 V.SetSize(vdim);
375 if (Function)
376 {
377 Function(transip, V);
378 }
379 else
380 {
381 TDFunction(transip, GetTime(), V);
382 }
383 if (Q)
384 {
385 V *= Q->Eval(T, ip, GetTime());
386 }
387}
388
390 : VectorCoefficient(dim), Coeff(dim), ownCoeff(dim)
391{
392 for (int i = 0; i < dim; i++)
393 {
394 Coeff[i] = NULL;
395 ownCoeff[i] = true;
396 }
397}
398
400{
401 for (int i = 0; i < vdim; i++)
402 {
403 if (Coeff[i]) { Coeff[i]->SetTime(t); }
404 }
406}
407
409{
410 if (ownCoeff[i]) { delete Coeff[i]; }
411 Coeff[i] = c;
412 ownCoeff[i] = own;
413}
414
416{
417 for (int i = 0; i < vdim; i++)
418 {
419 if (ownCoeff[i]) { delete Coeff[i]; }
420 }
421}
422
424 const IntegrationPoint &ip)
425{
426 V.SetSize(vdim);
427 for (int i = 0; i < vdim; i++)
428 {
429 V(i) = this->Eval(i, T, ip);
430 }
431}
432
434 const GridFunction *gf)
435 : VectorCoefficient ((gf) ? gf -> VectorDim() : 0)
436{
437 GridFunc = gf;
438}
439
441{
442 GridFunc = gf; vdim = (gf) ? gf -> VectorDim() : 0;
443}
444
446 const IntegrationPoint &ip)
447{
448 Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
449 if (T.mesh->GetNE() == gf_mesh->GetNE())
450 {
451 GridFunc->GetVectorValue(T, ip, V);
452 }
453 else
454 {
455 IntegrationPoint coarse_ip;
456 ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
457 GridFunc->GetVectorValue(*coarse_T, coarse_ip, V);
458 }
459}
460
463{
464 if (T.mesh == GridFunc->FESpace()->GetMesh())
465 {
466 GridFunc->GetVectorValues(T, ir, M);
467 }
468 else
469 {
470 VectorCoefficient::Eval(M, T, ir);
471 }
472}
473
478
480 const GridFunction *gf)
481 : VectorCoefficient((gf) ?
482 gf -> FESpace() -> GetMesh() -> SpaceDimension() : 0)
483{
484 GridFunc = gf;
485}
486
488{
489 GridFunc = gf; vdim = (gf) ?
490 gf -> FESpace() -> GetMesh() -> SpaceDimension() : 0;
491}
492
494 const IntegrationPoint &ip)
495{
496 Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
497 if (T.mesh->GetNE() == gf_mesh->GetNE())
498 {
499 GridFunc->GetGradient(T, V);
500 }
501 else
502 {
503 IntegrationPoint coarse_ip;
504 ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
505 GridFunc->GetGradient(*coarse_T, V);
506 }
507}
508
511{
512 if (T.mesh == GridFunc->FESpace()->GetMesh())
513 {
514 GridFunc->GetGradients(T, ir, M);
515 }
516 else
517 {
518 VectorCoefficient::Eval(M, T, ir);
519 }
520}
521
528
530{
531 GridFunc = gf; vdim = (gf) ? gf -> CurlDim() : 0;
532}
533
535 const IntegrationPoint &ip)
536{
537 Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
538 if (T.mesh->GetNE() == gf_mesh->GetNE())
539 {
540 GridFunc->GetCurl(T, V);
541 }
542 else
543 {
544 IntegrationPoint coarse_ip;
545 ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
546 GridFunc->GetCurl(*coarse_T, V);
547 }
548}
549
555
557 const IntegrationPoint &ip)
558{
559 Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
560 if (T.mesh->GetNE() == gf_mesh->GetNE())
561 {
562 return GridFunc->GetDivergence(T);
563 }
564 else
565 {
566 IntegrationPoint coarse_ip;
567 ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
568 return GridFunc->GetDivergence(*coarse_T);
569 }
570}
571
577
579{
580 dir = d_;
581 (*this).vdim = dir.Size();
582}
583
586{
587 V = dir;
588 d.SetTime(GetTime());
589 V *= d.EvalDelta(T, ip);
590}
591
593{
594 if (c) { c->SetTime(t); }
596}
597
599 const IntegrationPoint &ip)
600{
601 V.SetSize(vdim);
602 if (active_attr[T.Attribute-1])
603 {
604 c->SetTime(GetTime());
605 c->Eval(V, T, ip);
606 }
607 else
608 {
609 V = 0.0;
610 }
611}
612
615{
616 if (active_attr[T.Attribute-1])
617 {
618 c->SetTime(GetTime());
619 c->Eval(M, T, ir);
620 }
621 else
622 {
623 M.SetSize(vdim, ir.GetNPoints());
624 M = 0.0;
625 }
626}
627
629{
630 MFEM_VERIFY(qf.GetVDim() == height*width, "Wrong sizes.");
631 QuadratureSpaceBase &qspace = *qf.GetSpace();
632 const int ne = qspace.GetNE();
633 DenseMatrix values, matrix;
634 for (int iel = 0; iel < ne; ++iel)
635 {
636 qf.GetValues(iel, values);
637 const IntegrationRule &ir = qspace.GetIntRule(iel);
639 for (int iq = 0; iq < ir.Size(); ++iq)
640 {
641 const IntegrationPoint &ip = ir[iq];
642 T.SetIntPoint(&ip);
643 const int iq_p = qspace.GetPermutedIndex(iel, iq);
644 matrix.UseExternalData(&values(0, iq_p), height, width);
645 Eval(matrix, T, ip);
646 if (transpose) { matrix.Transpose(); }
647 }
648 }
649}
650
651void PWMatrixCoefficient::InitMap(const Array<int> & attr,
652 const Array<MatrixCoefficient*> & coefs)
653{
654 MFEM_VERIFY(attr.Size() == coefs.Size(),
655 "PWMatrixCoefficient: "
656 "Attribute and coefficient arrays have incompatible "
657 "dimensions.");
658
659 for (int i=0; i<attr.Size(); i++)
660 {
661 if (coefs[i] != NULL)
662 {
663 UpdateCoefficient(attr[i], *coefs[i]);
664 }
665 }
666}
667
669{
670 MFEM_VERIFY(coef.GetHeight() == height,
671 "PWMatrixCoefficient::UpdateCoefficient: "
672 "MatrixCoefficient has incompatible height.");
673 MFEM_VERIFY(coef.GetWidth() == width,
674 "PWMatrixCoefficient::UpdateCoefficient: "
675 "MatrixCoefficient has incompatible width.");
676 if (symmetric)
677 {
678 MFEM_VERIFY(coef.IsSymmetric(),
679 "PWMatrixCoefficient::UpdateCoefficient: "
680 "MatrixCoefficient has incompatible symmetry.");
681 }
682 pieces[attr] = &coef;
683}
684
686{
688
689 std::map<int, MatrixCoefficient*>::iterator p = pieces.begin();
690 for (; p != pieces.end(); p++)
691 {
692 if (p->second != NULL)
693 {
694 p->second->SetTime(t);
695 }
696 }
697}
698
700 const IntegrationPoint &ip)
701{
702 const int att = T.Attribute;
703 std::map<int, MatrixCoefficient*>::const_iterator p = pieces.find(att);
704 if (p != pieces.end())
705 {
706 if ( p->second != NULL)
707 {
708 p->second->Eval(K, T, ip);
709 return;
710 }
711 }
712
713 K.SetSize(height, width);
714 K = 0.0;
715}
716
718{
719 if (Q) { Q->SetTime(t); }
721}
722
724 const IntegrationPoint &ip)
725{
726 real_t x[3];
727 Vector transip(x, 3);
728
729 T.Transform(ip, transip);
730
731 K.SetSize(height, width);
732
733 if (symmetric) // Use SymmFunction (deprecated version)
734 {
735 MFEM_VERIFY(height == width && SymmFunction,
736 "MatrixFunctionCoefficient is not symmetric");
737
738 Vector Ksym((width * (width + 1)) / 2); // 1x1: 1, 2x2: 3, 3x3: 6
739
740 SymmFunction(transip, Ksym);
741
742 // Copy upper triangular values from Ksym to the full matrix K
743 int os = 0;
744 for (int i=0; i<height; ++i)
745 {
746 for (int j=i; j<width; ++j)
747 {
748 const real_t Kij = Ksym[j - i + os];
749 K(i,j) = Kij;
750 if (j != i) { K(j,i) = Kij; }
751 }
752
753 os += width - i;
754 }
755 }
756 else
757 {
758 if (Function)
759 {
760 Function(transip, K);
761 }
762 else if (TDFunction)
763 {
764 TDFunction(transip, GetTime(), K);
765 }
766 else
767 {
768 K = mat;
769 }
770 }
771
772 if (Q)
773 {
774 K *= Q->Eval(T, ip, GetTime());
775 }
776}
777
780 const IntegrationPoint &ip)
781{
782 MFEM_VERIFY(symmetric && height == width && SymmFunction,
783 "MatrixFunctionCoefficient is not symmetric");
784
785 real_t x[3];
786 Vector transip(x, 3);
787
788 T.Transform(ip, transip);
789
790 K.SetSize((width * (width + 1)) / 2); // 1x1: 1, 2x2: 3, 3x3: 6
791
792 if (SymmFunction)
793 {
794 SymmFunction(transip, K);
795 }
796
797 if (Q)
798 {
799 K *= Q->Eval(T, ip, GetTime());
800 }
801}
802
804{
805 const int vdim = qf.GetVDim();
806 MFEM_VERIFY(vdim == height*(height+1)/2, "Wrong sizes.");
807
808 QuadratureSpaceBase &qspace = *qf.GetSpace();
809 const int ne = qspace.GetNE();
810 DenseMatrix values;
812 for (int iel = 0; iel < ne; ++iel)
813 {
814 qf.GetValues(iel, values);
815 const IntegrationRule &ir = qspace.GetIntRule(iel);
817 for (int iq = 0; iq < ir.Size(); ++iq)
818 {
819 const IntegrationPoint &ip = ir[iq];
820 T.SetIntPoint(&ip);
821 matrix.UseExternalData(&values(0, iq), vdim);
822 Eval(matrix, T, ip);
823 }
824 }
825}
826
827
829 const IntegrationPoint &ip)
830{
832 Eval(mat, T, ip);
833 for (int j = 0; j < width; ++j)
834 {
835 for (int i = 0; i < height; ++ i)
836 {
837 K(i, j) = mat(i, j);
838 }
839 }
840}
841
847
850 const IntegrationPoint &ip)
851{
852 real_t x[3];
853 Vector transip(x, 3);
854
855 T.Transform(ip, transip);
856
857 K.SetSize(height);
858
859 if (Function)
860 {
861 Function(transip, K);
862 }
863 else if (TDFunction)
864 {
865 TDFunction(transip, GetTime(), K);
866 }
867 else
868 {
869 K = mat;
870 }
871
872 if (Q)
873 {
874 K *= Q->Eval(T, ip, GetTime());
875 }
876}
877
880{
881 Coeff.SetSize(height*width);
882 ownCoeff.SetSize(height*width);
883 for (int i = 0; i < (height*width); i++)
884 {
885 Coeff[i] = NULL;
886 ownCoeff[i] = true;
887 }
888}
889
891{
892 for (int i=0; i < height*width; i++)
893 {
894 if (Coeff[i]) { Coeff[i]->SetTime(t); }
895 }
897}
898
899void MatrixArrayCoefficient::Set(int i, int j, Coefficient * c, bool own)
900{
901 if (ownCoeff[i*width+j]) { delete Coeff[i*width+j]; }
902 Coeff[i*width+j] = c;
903 ownCoeff[i*width+j] = own;
904}
905
907{
908 for (int i=0; i < height*width; i++)
909 {
910 if (ownCoeff[i]) { delete Coeff[i]; }
911 }
912}
913
915 const IntegrationPoint &ip)
916{
917 K.SetSize(height, width);
918 for (int i = 0; i < height; i++)
919 {
920 for (int j = 0; j < width; j++)
921 {
922 K(i,j) = this->Eval(i, j, T, ip);
923 }
924 }
925}
926
928{
929 if (c) { c->SetTime(t); }
931}
932
934 const IntegrationPoint &ip)
935{
936 if (active_attr[T.Attribute-1])
937 {
938 c->SetTime(GetTime());
939 c->Eval(K, T, ip);
940 }
941 else
942 {
943 K.SetSize(height, width);
944 K = 0.0;
945 }
946}
947
949{
950 if (a) { a->SetTime(t); }
951 if (b) { b->SetTime(t); }
952 this->Coefficient::SetTime(t);
953}
954
956{
957 if (a) { a->SetTime(t); }
958 if (b) { b->SetTime(t); }
959 this->Coefficient::SetTime(t);
960}
961
963{
964 if (a) { a->SetTime(t); }
965 if (b) { b->SetTime(t); }
966 this->Coefficient::SetTime(t);
967}
968
970{
971 if (a) { a->SetTime(t); }
972 this->Coefficient::SetTime(t);
973}
974
977 : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
978{
979 MFEM_ASSERT(A.GetVDim() == B.GetVDim(),
980 "InnerProductCoefficient: "
981 "Arguments have incompatible dimensions.");
982}
983
985{
986 if (a) { a->SetTime(t); }
987 if (b) { b->SetTime(t); }
988 this->Coefficient::SetTime(t);
989}
990
992 const IntegrationPoint &ip)
993{
994 a->Eval(va, T, ip);
995 b->Eval(vb, T, ip);
996 return va * vb;
997}
998
1001 : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
1002{
1003 MFEM_ASSERT(A.GetVDim() == 2 && B.GetVDim() == 2,
1004 "VectorRotProductCoefficient: "
1005 "Arguments must have dimension equal to two.");
1006}
1007
1009{
1010 if (a) { a->SetTime(t); }
1011 if (b) { b->SetTime(t); }
1012 this->Coefficient::SetTime(t);
1013}
1014
1016 const IntegrationPoint &ip)
1017{
1018 a->Eval(va, T, ip);
1019 b->Eval(vb, T, ip);
1020 return va[0] * vb[1] - va[1] * vb[0];
1021}
1022
1024 : a(&A), ma(A.GetHeight(), A.GetWidth())
1025{
1026 MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1027 "DeterminantCoefficient: "
1028 "Argument must be a square matrix.");
1029}
1030
1032{
1033 if (a) { a->SetTime(t); }
1034 this->Coefficient::SetTime(t);
1035}
1036
1038 const IntegrationPoint &ip)
1039{
1040 a->Eval(ma, T, ip);
1041 return ma.Det();
1042}
1043
1046 ACoef(NULL), BCoef(NULL),
1047 A(dim), B(dim),
1048 alphaCoef(NULL), betaCoef(NULL),
1049 alpha(1.0), beta(1.0)
1050{
1051 A = 0.0; B = 0.0;
1052}
1053
1056 real_t alpha_, real_t beta_)
1057 : VectorCoefficient(A_.GetVDim()),
1058 ACoef(&A_), BCoef(&B_),
1059 A(A_.GetVDim()), B(A_.GetVDim()),
1060 alphaCoef(NULL), betaCoef(NULL),
1061 alpha(alpha_), beta(beta_)
1062{
1063 MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
1064 "VectorSumCoefficient: "
1065 "Arguments must have the same dimension.");
1066}
1067
1070 Coefficient &alpha_,
1071 Coefficient &beta_)
1072 : VectorCoefficient(A_.GetVDim()),
1073 ACoef(&A_), BCoef(&B_),
1074 A(A_.GetVDim()),
1075 B(A_.GetVDim()),
1076 alphaCoef(&alpha_),
1077 betaCoef(&beta_),
1078 alpha(0.0), beta(0.0)
1079{
1080 MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
1081 "VectorSumCoefficient: "
1082 "Arguments must have the same dimension.");
1083}
1084
1086{
1087 if (ACoef) { ACoef->SetTime(t); }
1088 if (BCoef) { BCoef->SetTime(t); }
1089 if (alphaCoef) { alphaCoef->SetTime(t); }
1090 if (betaCoef) { betaCoef->SetTime(t); }
1092}
1093
1095 const IntegrationPoint &ip)
1096{
1097 V.SetSize(A.Size());
1098 if ( ACoef) { ACoef->Eval(A, T, ip); }
1099 if ( BCoef) { BCoef->Eval(B, T, ip); }
1100 if (alphaCoef) { alpha = alphaCoef->Eval(T, ip); }
1101 if ( betaCoef) { beta = betaCoef->Eval(T, ip); }
1102 add(alpha, A, beta, B, V);
1103}
1104
1106 real_t A,
1108 : VectorCoefficient(B.GetVDim()), aConst(A), a(NULL), b(&B)
1109{}
1110
1116
1118{
1119 if (a) { a->SetTime(t); }
1120 if (b) { b->SetTime(t); }
1122}
1123
1125 const IntegrationPoint &ip)
1126{
1127 real_t sa = (a == NULL) ? aConst : a->Eval(T, ip);
1128 b->Eval(V, T, ip);
1129 V *= sa;
1130}
1131
1136
1138{
1139 if (a) { a->SetTime(t); }
1141}
1142
1144 const IntegrationPoint &ip)
1145{
1146 a->Eval(V, T, ip);
1147 real_t nv = V.Norml2();
1148 V *= (nv > tol) ? (1.0/nv) : 0.0;
1149}
1150
1154 : VectorCoefficient(3), a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
1155{
1156 MFEM_ASSERT(A.GetVDim() == 3 && B.GetVDim() == 3,
1157 "VectorCrossProductCoefficient: "
1158 "Arguments must have dimension equal to three.");
1159}
1160
1162{
1163 if (a) { a->SetTime(t); }
1164 if (b) { b->SetTime(t); }
1166}
1167
1169 const IntegrationPoint &ip)
1170{
1171 a->Eval(va, T, ip);
1172 b->Eval(vb, T, ip);
1173 V.SetSize(3);
1174 V[0] = va[1] * vb[2] - va[2] * vb[1];
1175 V[1] = va[2] * vb[0] - va[0] * vb[2];
1176 V[2] = va[0] * vb[1] - va[1] * vb[0];
1177}
1178
1181 : VectorCoefficient(A.GetHeight()), a(&A), b(&B),
1182 ma(A.GetHeight(), A.GetWidth()), vb(B.GetVDim())
1183{
1184 MFEM_ASSERT(A.GetWidth() == B.GetVDim(),
1185 "MatrixVectorProductCoefficient: "
1186 "Arguments have incompatible dimensions.");
1187}
1188
1190{
1191 if (a) { a->SetTime(t); }
1192 if (b) { b->SetTime(t); }
1194}
1195
1197 const IntegrationPoint &ip)
1198{
1199 a->Eval(ma, T, ip);
1200 b->Eval(vb, T, ip);
1201 V.SetSize(vdim);
1202 ma.Mult(vb, V);
1203}
1204
1206 const IntegrationPoint &ip)
1207{
1208 M.SetSize(dim);
1209 M = 0.0;
1210 for (int d=0; d<dim; d++) { M(d,d) = 1.0; }
1211}
1212
1215 real_t alpha_, real_t beta_)
1216 : MatrixCoefficient(A.GetHeight(), A.GetWidth()),
1217 a(&A), b(&B), alpha(alpha_), beta(beta_),
1218 ma(A.GetHeight(), A.GetWidth())
1219{
1220 MFEM_ASSERT(A.GetHeight() == B.GetHeight() && A.GetWidth() == B.GetWidth(),
1221 "MatrixSumCoefficient: "
1222 "Arguments must have the same dimensions.");
1223}
1224
1226{
1227 if (a) { a->SetTime(t); }
1228 if (b) { b->SetTime(t); }
1230}
1231
1233 const IntegrationPoint &ip)
1234{
1235 b->Eval(M, T, ip);
1236 if ( beta != 1.0 ) { M *= beta; }
1237 a->Eval(ma, T, ip);
1238 M.Add(alpha, ma);
1239}
1240
1243 : MatrixCoefficient(A.GetHeight(), B.GetWidth()),
1244 a(&A), b(&B),
1245 ma(A.GetHeight(), A.GetWidth()),
1246 mb(B.GetHeight(), B.GetWidth())
1247{
1248 MFEM_ASSERT(A.GetWidth() == B.GetHeight(),
1249 "MatrixProductCoefficient: "
1250 "Arguments must have compatible dimensions.");
1251}
1252
1254 const IntegrationPoint &ip)
1255{
1256 a->Eval(ma, T, ip);
1257 b->Eval(mb, T, ip);
1258 Mult(ma, mb, M);
1259}
1260
1262 real_t A,
1264 : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(A), a(NULL), b(&B)
1265{}
1266
1268 Coefficient &A,
1270 : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(0.0), a(&A), b(&B)
1271{}
1272
1274{
1275 if (a) { a->SetTime(t); }
1276 if (b) { b->SetTime(t); }
1278}
1279
1282 const IntegrationPoint &ip)
1283{
1284 real_t sa = (a == NULL) ? aConst : a->Eval(T, ip);
1285 b->Eval(M, T, ip);
1286 M *= sa;
1287}
1288
1292
1294{
1295 if (a) { a->SetTime(t); }
1297}
1298
1301 const IntegrationPoint &ip)
1302{
1303 a->Eval(M, T, ip);
1304 M.Transpose();
1305}
1306
1308 : MatrixCoefficient(A.GetHeight(), A.GetWidth()), a(&A)
1309{
1310 MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1311 "InverseMatrixCoefficient: "
1312 "Argument must be a square matrix.");
1313}
1314
1316{
1317 if (a) { a->SetTime(t); }
1319}
1320
1323 const IntegrationPoint &ip)
1324{
1325 a->Eval(M, T, ip);
1326 M.Invert();
1327}
1328
1331 : MatrixCoefficient(A.GetVDim(), B.GetVDim()), a(&A), b(&B),
1332 va(A.GetVDim()), vb(B.GetVDim())
1333{}
1334
1336{
1337 if (a) { a->SetTime(t); }
1338 if (b) { b->SetTime(t); }
1340}
1341
1343 const IntegrationPoint &ip)
1344{
1345 a->Eval(va, T, ip);
1346 b->Eval(vb, T, ip);
1347 M.SetSize(va.Size(), vb.Size());
1348 for (int i=0; i<va.Size(); i++)
1349 {
1350 for (int j=0; j<vb.Size(); j++)
1351 {
1352 M(i, j) = va[i] * vb[j];
1353 }
1354 }
1355}
1356
1358 : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(A), a(NULL), k(&K),
1359 vk(K.GetVDim())
1360{}
1361
1364 : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(0.0), a(&A), k(&K),
1365 vk(K.GetVDim())
1366{}
1367
1369{
1370 if (a) { a->SetTime(t); }
1371 if (k) { k->SetTime(t); }
1373}
1374
1376 const IntegrationPoint &ip)
1377{
1378 k->Eval(vk, T, ip);
1379 M.SetSize(vk.Size(), vk.Size());
1380 M = 0.0;
1381 real_t k2 = vk*vk;
1382 for (int i=0; i<vk.Size(); i++)
1383 {
1384 M(i, i) = k2;
1385 for (int j=0; j<vk.Size(); j++)
1386 {
1387 M(i, j) -= vk[i] * vk[j];
1388 }
1389 }
1390 M *= ((a == NULL ) ? aConst : a->Eval(T, ip) );
1391}
1392
1394 const IntegrationRule *irs[])
1395{
1396 real_t norm = 0.0;
1398
1399 for (int i = 0; i < mesh.GetNE(); i++)
1400 {
1401 tr = mesh.GetElementTransformation(i);
1402 const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1403 for (int j = 0; j < ir.GetNPoints(); j++)
1404 {
1405 const IntegrationPoint &ip = ir.IntPoint(j);
1406 tr->SetIntPoint(&ip);
1407 real_t val = fabs(coeff.Eval(*tr, ip));
1408 if (p < infinity())
1409 {
1410 norm += ip.weight * tr->Weight() * pow(val, p);
1411 }
1412 else
1413 {
1414 if (norm < val)
1415 {
1416 norm = val;
1417 }
1418 }
1419 }
1420 }
1421 return norm;
1422}
1423
1425 const IntegrationRule *irs[])
1426{
1427 real_t norm = 0.0;
1429 int vdim = coeff.GetVDim();
1430 Vector vval(vdim);
1431 real_t val;
1432
1433 for (int i = 0; i < mesh.GetNE(); i++)
1434 {
1435 tr = mesh.GetElementTransformation(i);
1436 const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1437 for (int j = 0; j < ir.GetNPoints(); j++)
1438 {
1439 const IntegrationPoint &ip = ir.IntPoint(j);
1440 tr->SetIntPoint(&ip);
1441 coeff.Eval(vval, *tr, ip);
1442 if (p < infinity())
1443 {
1444 for (int idim(0); idim < vdim; ++idim)
1445 {
1446 norm += ip.weight * tr->Weight() * pow(fabs( vval(idim) ), p);
1447 }
1448 }
1449 else
1450 {
1451 for (int idim(0); idim < vdim; ++idim)
1452 {
1453 val = fabs(vval(idim));
1454 if (norm < val)
1455 {
1456 norm = val;
1457 }
1458 }
1459 }
1460 }
1461 }
1462
1463 return norm;
1465
1467 const IntegrationRule *irs[])
1468{
1469 real_t norm = LpNormLoop(p, coeff, mesh, irs);
1470
1471 if (p < infinity())
1472 {
1473 // negative quadrature weights may cause norm to be negative
1474 if (norm < 0.0)
1475 {
1476 norm = -pow(-norm, 1.0/p);
1477 }
1478 else
1479 {
1480 norm = pow(norm, 1.0/p);
1481 }
1482 }
1483
1484 return norm;
1485}
1486
1488 const IntegrationRule *irs[])
1489{
1490 real_t norm = LpNormLoop(p, coeff, mesh, irs);
1491
1492 if (p < infinity())
1493 {
1494 // negative quadrature weights may cause norm to be negative
1495 if (norm < 0.0)
1496 {
1497 norm = -pow(-norm, 1.0/p);
1498 }
1499 else
1500 {
1501 norm = pow(norm, 1.0/p);
1502 }
1503 }
1504
1505 return norm;
1506}
1507
1508#ifdef MFEM_USE_MPI
1510 const IntegrationRule *irs[])
1511{
1512 real_t loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1513 real_t glob_norm = 0;
1514
1515 MPI_Comm comm = pmesh.GetComm();
1516
1517 if (p < infinity())
1518 {
1519 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
1520 comm);
1521
1522 // negative quadrature weights may cause norm to be negative
1523 if (glob_norm < 0.0)
1524 {
1525 glob_norm = -pow(-glob_norm, 1.0/p);
1526 }
1527 else
1528 {
1529 glob_norm = pow(glob_norm, 1.0/p);
1530 }
1531 }
1532 else
1533 {
1534 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
1535 comm);
1536 }
1537
1538 return glob_norm;
1539}
1540
1542 const IntegrationRule *irs[])
1543{
1544 real_t loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1545 real_t glob_norm = 0;
1546
1547 MPI_Comm comm = pmesh.GetComm();
1548
1549 if (p < infinity())
1550 {
1551 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
1552 comm);
1553
1554 // negative quadrature weights may cause norm to be negative
1555 if (glob_norm < 0.0)
1556 {
1557 glob_norm = -pow(-glob_norm, 1.0/p);
1558 }
1559 else
1560 {
1561 glob_norm = pow(glob_norm, 1.0/p);
1562 }
1563 }
1564 else
1565 {
1566 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
1567 comm);
1568 }
1569
1570 return glob_norm;
1571}
1572#endif
1573
1577
1579{
1580 MFEM_VERIFY(index_ >= 0, "Index must be >= 0");
1581 MFEM_VERIFY(index_ < QuadF.GetVDim(),
1582 "Index must be < QuadratureFunction length");
1583 index = index_;
1584
1585 MFEM_VERIFY(length_ > 0, "Length must be > 0");
1586 MFEM_VERIFY(length_ <= QuadF.GetVDim() - index,
1587 "Length must be <= (QuadratureFunction length - index)");
1588
1589 vdim = length_;
1590}
1591
1594 const IntegrationPoint &ip)
1595{
1596 QuadF.HostRead();
1597
1598 const int el_idx = QuadF.GetSpace()->GetEntityIndex(T);
1599 // Handle the case of "interior boundary elements" and FaceQuadratureSpace
1600 // with FaceType::Boundary.
1601 if (el_idx < 0) { V = 0.0; return; }
1602
1603 const int ip_idx = QuadF.GetSpace()->GetPermutedIndex(el_idx, ip.index);
1604
1605 if (index == 0 && vdim == QuadF.GetVDim())
1606 {
1607 QuadF.GetValues(el_idx, ip_idx, V);
1608 }
1609 else
1610 {
1611 Vector temp;
1612 QuadF.GetValues(el_idx, ip_idx, temp);
1613 V.SetSize(vdim);
1614 for (int i = 0; i < vdim; i++)
1615 {
1616 V(i) = temp(index + i);
1617 }
1618 }
1619
1620 return;
1621}
1622
1627
1629 const QuadratureFunction &qf) : QuadF(qf)
1630{
1631 MFEM_VERIFY(qf.GetVDim() == 1, "QuadratureFunction's vdim must be 1");
1632}
1633
1635 const IntegrationPoint &ip)
1636{
1637 QuadF.HostRead();
1638 Vector temp(1);
1639 const int el_idx = QuadF.GetSpace()->GetEntityIndex(T);
1640 // Handle the case of "interior boundary elements" and FaceQuadratureSpace
1641 // with FaceType::Boundary.
1642 if (el_idx < 0) { return 0.0; }
1643 const int ip_idx = QuadF.GetSpace()->GetPermutedIndex(el_idx, ip.index);
1644 QuadF.GetValues(el_idx, ip_idx, temp);
1645 return temp[0];
1646}
1647
1649{
1650 qf = QuadF;
1651}
1652
1653
1656 : Vector(), storage(storage_), vdim(0), qs(qs_), qf(NULL)
1657{
1658 UseDevice(true);
1659}
1660
1663 CoefficientStorage storage_)
1664 : CoefficientVector(qs_, storage_)
1665{
1666 if (coeff == NULL)
1667 {
1668 SetConstant(1.0);
1669 }
1670 else
1671 {
1672 Project(*coeff);
1673 }
1674}
1675
1678 CoefficientStorage storage_)
1679 : CoefficientVector(qs_, storage_)
1680{
1681 Project(coeff);
1682}
1683
1686 CoefficientStorage storage_)
1687 : CoefficientVector(qs_, storage_)
1688{
1689 Project(coeff);
1690}
1691
1694 CoefficientStorage storage_)
1695 : CoefficientVector(qs_, storage_)
1696{
1697 Project(coeff);
1698}
1699
1701{
1702 vdim = 1;
1703 if (auto *const_coeff = dynamic_cast<ConstantCoefficient*>(&coeff))
1704 {
1705 SetConstant(const_coeff->constant);
1706 }
1707 else if (auto *qf_coeff = dynamic_cast<QuadratureFunctionCoefficient*>(&coeff))
1708 {
1709 MakeRef(qf_coeff->GetQuadFunction());
1710 }
1711 else
1712 {
1713 if (qf == nullptr) { qf = new QuadratureFunction(qs); }
1714 qf->SetVDim(1);
1715 coeff.Project(*qf);
1716 Vector::MakeRef(*qf, 0, qf->Size());
1717 }
1718}
1719
1721{
1722 vdim = coeff.GetVDim();
1723 if (auto *const_coeff = dynamic_cast<VectorConstantCoefficient*>(&coeff))
1724 {
1725 SetConstant(const_coeff->GetVec());
1726 }
1727 else if (auto *qf_coeff =
1728 dynamic_cast<VectorQuadratureFunctionCoefficient*>(&coeff))
1729 {
1730 MakeRef(qf_coeff->GetQuadFunction());
1731 }
1732 else
1733 {
1734 if (qf == nullptr) { qf = new QuadratureFunction(qs, vdim); }
1735 qf->SetVDim(vdim);
1736 coeff.Project(*qf);
1737 Vector::MakeRef(*qf, 0, qf->Size());
1738 }
1739}
1740
1742{
1743 if (auto *const_coeff = dynamic_cast<MatrixConstantCoefficient*>(&coeff))
1744 {
1745 SetConstant(const_coeff->GetMatrix());
1746 }
1747 else if (auto *const_sym_coeff =
1748 dynamic_cast<SymmetricMatrixConstantCoefficient*>(&coeff))
1749 {
1750 SetConstant(const_sym_coeff->GetMatrix());
1751 }
1752 else
1753 {
1754 auto *sym_coeff = dynamic_cast<SymmetricMatrixCoefficient*>(&coeff);
1755 const bool sym = sym_coeff && (storage & CoefficientStorage::SYMMETRIC);
1756 const int height = coeff.GetHeight();
1757 const int width = coeff.GetWidth();
1758 vdim = sym ? height*(height + 1)/2 : width*height;
1759
1760 if (qf == nullptr) { qf = new QuadratureFunction(qs, vdim); }
1761 qf->SetVDim(vdim);
1762 if (sym) { sym_coeff->ProjectSymmetric(*qf); }
1763 else { coeff.Project(*qf, transpose); }
1764 Vector::MakeRef(*qf, 0, qf->Size());
1765 }
1766}
1767
1769{
1770 Project(coeff, true);
1771}
1772
1774{
1775 vdim = qf_.GetVDim();
1776 const QuadratureSpaceBase *qs2 = qf_.GetSpace();
1777 MFEM_CONTRACT_VAR(qs2); // qs2 used only for asserts
1778 MFEM_VERIFY(qs2 != NULL, "Invalid QuadratureSpace.")
1779 MFEM_VERIFY(qs2->GetMesh() == qs.GetMesh(), "Meshes differ.");
1780 MFEM_VERIFY(qs2->GetOrder() == qs.GetOrder(), "Orders differ.");
1781 Vector::MakeRef(const_cast<QuadratureFunction&>(qf_), 0, qf_.Size());
1782}
1783
1785{
1786 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1787 vdim = 1;
1788 SetSize(nq);
1789 Vector::operator=(constant);
1790}
1791
1793{
1794 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1795 vdim = constant.Size();
1796 SetSize(nq*vdim);
1797 for (int iq = 0; iq < nq; ++iq)
1798 {
1799 for (int vd = 0; vd<vdim; ++vd)
1800 {
1801 (*this)[vd + iq*vdim] = constant[vd];
1802 }
1803 }
1804}
1805
1807{
1808 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1809 const int width = constant.Width();
1810 const int height = constant.Height();
1811 vdim = width*height;
1812 SetSize(nq*vdim);
1813 for (int iq = 0; iq < nq; ++iq)
1814 {
1815 for (int j = 0; j < width; ++j)
1816 {
1817 for (int i = 0; i < height; ++i)
1818 {
1819 (*this)[i + j*height + iq*vdim] = constant(i, j);
1820 }
1821 }
1822 }
1823}
1824
1826{
1827 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1828 const int height = constant.Height();
1829 const bool sym = storage & CoefficientStorage::SYMMETRIC;
1830 vdim = sym ? height*(height + 1)/2 : height*height;
1831 SetSize(nq*vdim);
1832 for (int iq = 0; iq < nq; ++iq)
1833 {
1834 for (int vd = 0; vd < vdim; ++vd)
1835 {
1836 const real_t value = sym ? constant.GetData()[vd] : constant(vd % height,
1837 vd / height);
1838 (*this)[vd + iq*vdim] = value;
1839 }
1840 }
1841}
1842
1843int CoefficientVector::GetVDim() const { return vdim; }
1844
1849
1850}
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Class to represent a coefficient evaluated at quadrature points.
void SetConstant(real_t constant)
Set this vector to the given constant.
int vdim
Number of values per quadrature point.
void Project(Coefficient &coeff)
Evaluate the given Coefficient at the quadrature points defined by qs.
int GetVDim() const
Return the number of values per quadrature point.
CoefficientVector(QuadratureSpaceBase &qs_, CoefficientStorage storage_=CoefficientStorage::FULL)
Create an empty CoefficientVector.
QuadratureFunction * qf
Internal QuadratureFunction (owned, may be NULL).
CoefficientStorage storage
Storage optimizations (see CoefficientStorage).
QuadratureSpaceBase & qs
Associated QuadratureSpaceBase.
void ProjectTranspose(MatrixCoefficient &coeff)
Project the transpose of coeff.
void MakeRef(const QuadratureFunction &qf_)
Make this vector a reference to the given QuadratureFunction.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
real_t GetTime()
Get the time for time dependent coefficients.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
A coefficient that is constant across space and time.
void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf with the constant value.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
CrossCrossCoefficient(real_t A, VectorCoefficient &K)
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
CurlGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector curl coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
real_t Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
void SetDeltaCenter(const Vector &center)
Set the center location of the delta function.
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
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:220
void Transpose()
(*this) = (*this)^t
void GetColumnReference(int c, Vector &col)
Definition densemat.hpp:312
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
void Invert()
Replaces the current matrix with its inverse.
Definition densemat.cpp:726
void UseExternalData(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition densemat.hpp:80
void Add(const real_t c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition densemat.cpp:628
real_t Det() const
Definition densemat.cpp:535
void SetSize(int s)
Change the size of the DenseSymmetricMatrix to s x s.
Definition symmat.cpp:32
void UseExternalData(real_t *d, int s)
Change the data array and the size of the DenseSymmetricMatrix.
Definition symmat.hpp:48
real_t * GetData() const
Returns the matrix data array.
Definition symmat.hpp:80
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the determinant coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
DivergenceGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the scalar divergence coefficient at ip.
const Mesh * mesh
The Mesh object containing the element.
Definition eltrans.hpp:84
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:162
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:93
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
std::function< real_t(const Vector &)> Function
std::function< real_t(const Vector &, real_t)> TDFunction
void SetGridFunction(const GridFunction *gf)
Set the scalar grid function.
GradientGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a scalar grid function gf. The grid function is not owned by the coeff...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the gradient vector coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition gridfunc.cpp:449
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition gridfunc.cpp:476
real_t GetDivergence(ElementTransformation &tr) const
void GetCurl(ElementTransformation &tr, Vector &curl) const
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition gridfunc.cpp:714
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
void SetTime(real_t t)
Set the time for internally stored coefficients.
IsoparametricTransformation Transf
Definition eltrans.hpp:467
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition eltrans.cpp:543
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
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
A standard isoparametric element transformation.
Definition eltrans.hpp:363
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition eltrans.hpp:402
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition eltrans.cpp:379
void Set(int i, int j, Coefficient *c, bool own=true)
Set the coefficient located at (i,j) in the matrix. By default by default this will take ownership of...
MatrixArrayCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
real_t Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
void SetTime(real_t t)
Set the time for internally stored coefficients.
Base class for Matrix Coefficients that optionally depend on time and space.
virtual void Project(QuadratureFunction &qf, bool transpose=false)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points....
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t GetTime()
Get the time for time dependent coefficients.
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 GetWidth() const
Get the width of the matrix.
int GetHeight() const
Get the height of the matrix.
A matrix coefficient that is constant in space and time.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
MatrixSumCoefficient(MatrixCoefficient &A, MatrixCoefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Mesh data type.
Definition mesh.hpp:56
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition mesh.cpp:7316
const CoarseFineTransformations & GetRefinementTransforms() const
Definition mesh.cpp:11082
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
NormalizedVectorCoefficient(VectorCoefficient &A, real_t tol=1e-6)
Return a vector normalized to a length of one.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void UpdateCoefficient(int attr, Coefficient &coef)
Replace a single Coefficient for a particular attribute.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
QuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition qfunction.hpp:78
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition qfunction.hpp:82
void ProjectGridFunction(const GridFunction &gf)
Evaluate a grid function at each quadrature point.
void GetValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
int GetVDim() const
Get the vector dimension.
Definition qfunction.hpp:75
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition qspace.hpp:26
int GetNE() const
Return the number of entities.
Definition qspace.hpp:69
virtual ElementTransformation * GetTransformation(int idx)=0
Get the (element or face) transformation of entity idx.
int GetSize() const
Return the total number of quadrature points.
Definition qspace.hpp:63
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition qspace.hpp:81
virtual int GetEntityIndex(const ElementTransformation &T) const =0
Returns the index in the quadrature space of the entity associated with the transformation T.
virtual int GetPermutedIndex(int idx, int iq) const =0
Returns the permuted index of the iq quadrature point in entity idx.
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition qspace.hpp:66
Mesh * GetMesh() const
Returns the mesh.
Definition qspace.hpp:72
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetTime(real_t t)
Set the time for internally stored coefficients.
ScalarMatrixProductCoefficient(real_t A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
ScalarVectorProductCoefficient(real_t A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Base class for symmetric matrix coefficients that optionally depend on time and space.
DenseSymmetricMatrix mat
Internal matrix used when evaluating this coefficient as a DenseMatrix.
virtual void ProjectSymmetric(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result as ...
A matrix coefficient that is constant in space and time.
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set().
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
real_t Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t GetTime()
Get the time for time dependent coefficients.
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 Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Vector coefficient that is constant in space and time.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorCrossProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two coefficients. Result is A x B.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetDirection(const Vector &d_)
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 Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
VectorQuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
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 ...
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
void SetComponent(int index_, int length_)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Vector data type.
Definition vector.hpp:80
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:478
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:832
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:139
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
Vector & operator=(const real_t *v)
Copy Size() entries from v.
Definition vector.cpp:118
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:602
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
Mesh * GetMesh(int type)
Definition ex29.cpp:218
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:45
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:316
real_t ComputeLpNorm(real_t p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
CoefficientStorage
Flags that determine what storage optimizations to use in CoefficientVector.
@ SYMMETRIC
Store the triangular part of symmetric matrices.
@ CONSTANTS
Store constants using only vdim entries.
float real_t
Definition config.hpp:43
real_t LpNormLoop(real_t p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
ElementTransformation * RefinedToCoarse(Mesh &coarse_mesh, const ElementTransformation &T, const IntegrationPoint &ip, IntegrationPoint &coarse_ip)
real_t p(const Vector &x, real_t t)
RefCoord t[3]
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:72
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:74
DenseTensor point_matrices[Geometry::NumGeom]
Definition ncmesh.hpp:78
Helper struct to convert a C++ type to an MPI type.