MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
coefficient.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// 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 qf.HostWrite();
811 DenseMatrix values;
813 for (int iel = 0; iel < ne; ++iel)
814 {
815 qf.GetValues(iel, values);
816 const IntegrationRule &ir = qspace.GetIntRule(iel);
818 for (int iq = 0; iq < ir.Size(); ++iq)
819 {
820 const IntegrationPoint &ip = ir[iq];
821 T.SetIntPoint(&ip);
822 matrix.UseExternalData(&values(0, iq), height);
823 Eval(matrix, T, ip);
824 }
825 }
826}
827
828
830 const IntegrationPoint &ip)
831{
832 Eval(mat_aux, 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_aux(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
929{
930 Coeff.SetSize(height);
931 ownCoeff.SetSize(height);
932 for (int i = 0; i < height; i++)
933 {
934 Coeff[i] = NULL;
935 ownCoeff[i] = true;
936 }
937}
938
940{
941 for (int i=0; i < height; i++)
942 {
943 if (Coeff[i]) { Coeff[i]->SetTime(t); }
944 }
946}
947
949{
950 MFEM_ASSERT(i < height && i >= 0, "Row "
951 << i << " does not exist. " <<
952 "Matrix height = " << height << ".");
953 if (ownCoeff[i]) { delete Coeff[i]; }
954 Coeff[i] = c;
955 ownCoeff[i] = own;
956}
957
959{
960 for (int i=0; i < height; i++)
961 {
962 if (ownCoeff[i]) { delete Coeff[i]; }
963 }
964}
965
968 const IntegrationPoint &ip)
969{
970 MFEM_ASSERT(i < height && i >= 0, "Row "
971 << i << " does not exist. " <<
972 "Matrix height = " << height << ".");
973 if (Coeff[i])
974 {
975 Coeff[i] -> Eval(V, T, ip);
976 }
977 else
978 {
979 V = 0.0;
980 }
981}
982
985 const IntegrationPoint &ip)
986{
987 K.SetSize(height, width);
988 Vector V(width);
989 for (int i = 0; i < height; i++)
990 {
991 this->Eval(i, V, T, ip);
992 K.SetRow(i, V);
993 }
994}
995
997{
998 if (c) { c->SetTime(t); }
1000}
1001
1003 const IntegrationPoint &ip)
1004{
1005 if (active_attr[T.Attribute-1])
1006 {
1007 c->SetTime(GetTime());
1008 c->Eval(K, T, ip);
1009 }
1010 else
1011 {
1012 K.SetSize(height, width);
1013 K = 0.0;
1014 }
1015}
1016
1018{
1019 if (a) { a->SetTime(t); }
1020 if (b) { b->SetTime(t); }
1021 this->Coefficient::SetTime(t);
1022}
1023
1025{
1026 if (a) { a->SetTime(t); }
1027 if (b) { b->SetTime(t); }
1028 this->Coefficient::SetTime(t);
1029}
1030
1032{
1033 if (a) { a->SetTime(t); }
1034 if (b) { b->SetTime(t); }
1035 this->Coefficient::SetTime(t);
1036}
1037
1039{
1040 if (a) { a->SetTime(t); }
1041 this->Coefficient::SetTime(t);
1042}
1043
1046 : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
1047{
1048 MFEM_ASSERT(A.GetVDim() == B.GetVDim(),
1049 "InnerProductCoefficient: "
1050 "Arguments have incompatible dimensions.");
1051}
1052
1054{
1055 if (a) { a->SetTime(t); }
1056 if (b) { b->SetTime(t); }
1057 this->Coefficient::SetTime(t);
1058}
1059
1061 const IntegrationPoint &ip)
1062{
1063 a->Eval(va, T, ip);
1064 b->Eval(vb, T, ip);
1065 return va * vb;
1066}
1067
1070 : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
1071{
1072 MFEM_ASSERT(A.GetVDim() == 2 && B.GetVDim() == 2,
1073 "VectorRotProductCoefficient: "
1074 "Arguments must have dimension equal to two.");
1075}
1076
1078{
1079 if (a) { a->SetTime(t); }
1080 if (b) { b->SetTime(t); }
1081 this->Coefficient::SetTime(t);
1082}
1083
1085 const IntegrationPoint &ip)
1086{
1087 a->Eval(va, T, ip);
1088 b->Eval(vb, T, ip);
1089 return va[0] * vb[1] - va[1] * vb[0];
1090}
1091
1093 : a(&A), ma(A.GetHeight(), A.GetWidth())
1094{
1095 MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1096 "DeterminantCoefficient: "
1097 "Argument must be a square matrix.");
1098}
1099
1101{
1102 if (a) { a->SetTime(t); }
1103 this->Coefficient::SetTime(t);
1104}
1105
1107 const IntegrationPoint &ip)
1108{
1109 a->Eval(ma, T, ip);
1110 return ma.Det();
1111}
1112
1114 : a(&A), ma(A.GetHeight(), A.GetWidth())
1115{
1116 MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1117 "TraceCoefficient: "
1118 "Argument must be a square matrix.");
1119}
1120
1122{
1123 if (a) { a->SetTime(t); }
1124 this->Coefficient::SetTime(t);
1125}
1126
1128 const IntegrationPoint &ip)
1129{
1130 a->Eval(ma, T, ip);
1131 return ma.Trace();
1132}
1133
1136 ACoef(NULL), BCoef(NULL),
1137 A(dim), B(dim),
1138 alphaCoef(NULL), betaCoef(NULL),
1139 alpha(1.0), beta(1.0)
1140{
1141 A = 0.0; B = 0.0;
1142}
1143
1146 real_t alpha_, real_t beta_)
1147 : VectorCoefficient(A_.GetVDim()),
1148 ACoef(&A_), BCoef(&B_),
1149 A(A_.GetVDim()), B(A_.GetVDim()),
1150 alphaCoef(NULL), betaCoef(NULL),
1151 alpha(alpha_), beta(beta_)
1152{
1153 MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
1154 "VectorSumCoefficient: "
1155 "Arguments must have the same dimension.");
1156}
1157
1160 Coefficient &alpha_,
1161 Coefficient &beta_)
1162 : VectorCoefficient(A_.GetVDim()),
1163 ACoef(&A_), BCoef(&B_),
1164 A(A_.GetVDim()),
1165 B(A_.GetVDim()),
1166 alphaCoef(&alpha_),
1167 betaCoef(&beta_),
1168 alpha(0.0), beta(0.0)
1169{
1170 MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
1171 "VectorSumCoefficient: "
1172 "Arguments must have the same dimension.");
1173}
1174
1176{
1177 if (ACoef) { ACoef->SetTime(t); }
1178 if (BCoef) { BCoef->SetTime(t); }
1179 if (alphaCoef) { alphaCoef->SetTime(t); }
1180 if (betaCoef) { betaCoef->SetTime(t); }
1182}
1183
1185 const IntegrationPoint &ip)
1186{
1187 V.SetSize(A.Size());
1188 if ( ACoef) { ACoef->Eval(A, T, ip); }
1189 if ( BCoef) { BCoef->Eval(B, T, ip); }
1190 if (alphaCoef) { alpha = alphaCoef->Eval(T, ip); }
1191 if ( betaCoef) { beta = betaCoef->Eval(T, ip); }
1192 add(alpha, A, beta, B, V);
1193}
1194
1196 real_t A,
1198 : VectorCoefficient(B.GetVDim()), aConst(A), a(NULL), b(&B)
1199{}
1200
1206
1208{
1209 if (a) { a->SetTime(t); }
1210 if (b) { b->SetTime(t); }
1212}
1213
1215 const IntegrationPoint &ip)
1216{
1217 real_t sa = (a == NULL) ? aConst : a->Eval(T, ip);
1218 b->Eval(V, T, ip);
1219 V *= sa;
1220}
1221
1226
1228{
1229 if (a) { a->SetTime(t); }
1231}
1232
1234 const IntegrationPoint &ip)
1235{
1236 a->Eval(V, T, ip);
1237 real_t nv = V.Norml2();
1238 V *= (nv > tol) ? (1.0/nv) : 0.0;
1239}
1240
1244 : VectorCoefficient(3), a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
1245{
1246 MFEM_ASSERT(A.GetVDim() == 3 && B.GetVDim() == 3,
1247 "VectorCrossProductCoefficient: "
1248 "Arguments must have dimension equal to three.");
1249}
1250
1252{
1253 if (a) { a->SetTime(t); }
1254 if (b) { b->SetTime(t); }
1256}
1257
1259 const IntegrationPoint &ip)
1260{
1261 a->Eval(va, T, ip);
1262 b->Eval(vb, T, ip);
1263 V.SetSize(3);
1264 V[0] = va[1] * vb[2] - va[2] * vb[1];
1265 V[1] = va[2] * vb[0] - va[0] * vb[2];
1266 V[2] = va[0] * vb[1] - va[1] * vb[0];
1267}
1268
1271 : VectorCoefficient(A.GetHeight()), a(&A), b(&B),
1272 ma(A.GetHeight(), A.GetWidth()), vb(B.GetVDim())
1273{
1274 MFEM_ASSERT(A.GetWidth() == B.GetVDim(),
1275 "MatrixVectorProductCoefficient: "
1276 "Arguments have incompatible dimensions.");
1277}
1278
1280{
1281 if (a) { a->SetTime(t); }
1282 if (b) { b->SetTime(t); }
1284}
1285
1287 const IntegrationPoint &ip)
1288{
1289 a->Eval(ma, T, ip);
1290 b->Eval(vb, T, ip);
1291 V.SetSize(vdim);
1292 ma.Mult(vb, V);
1293}
1294
1296 const IntegrationPoint &ip)
1297{
1298 M.SetSize(dim);
1299 M = 0.0;
1300 for (int d=0; d<dim; d++) { M(d,d) = 1.0; }
1301}
1302
1305 real_t alpha_, real_t beta_)
1306 : MatrixCoefficient(A.GetHeight(), A.GetWidth()),
1307 a(&A), b(&B), alpha(alpha_), beta(beta_),
1308 ma(A.GetHeight(), A.GetWidth())
1309{
1310 MFEM_ASSERT(A.GetHeight() == B.GetHeight() && A.GetWidth() == B.GetWidth(),
1311 "MatrixSumCoefficient: "
1312 "Arguments must have the same dimensions.");
1313}
1314
1316{
1317 if (a) { a->SetTime(t); }
1318 if (b) { b->SetTime(t); }
1320}
1321
1323 const IntegrationPoint &ip)
1324{
1325 b->Eval(M, T, ip);
1326 if ( beta != 1.0 ) { M *= beta; }
1327 a->Eval(ma, T, ip);
1328 M.Add(alpha, ma);
1329}
1330
1333 : MatrixCoefficient(A.GetHeight(), B.GetWidth()),
1334 a(&A), b(&B),
1335 ma(A.GetHeight(), A.GetWidth()),
1336 mb(B.GetHeight(), B.GetWidth())
1337{
1338 MFEM_ASSERT(A.GetWidth() == B.GetHeight(),
1339 "MatrixProductCoefficient: "
1340 "Arguments must have compatible dimensions.");
1341}
1342
1344 const IntegrationPoint &ip)
1345{
1346 a->Eval(ma, T, ip);
1347 b->Eval(mb, T, ip);
1348 Mult(ma, mb, M);
1349}
1350
1352 real_t A,
1354 : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(A), a(NULL), b(&B)
1355{}
1356
1358 Coefficient &A,
1360 : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(0.0), a(&A), b(&B)
1361{}
1362
1364{
1365 if (a) { a->SetTime(t); }
1366 if (b) { b->SetTime(t); }
1368}
1369
1372 const IntegrationPoint &ip)
1373{
1374 real_t sa = (a == NULL) ? aConst : a->Eval(T, ip);
1375 b->Eval(M, T, ip);
1376 M *= sa;
1377}
1378
1382
1384{
1385 if (a) { a->SetTime(t); }
1387}
1388
1391 const IntegrationPoint &ip)
1392{
1393 a->Eval(M, T, ip);
1394 M.Transpose();
1395}
1396
1398 : MatrixCoefficient(A.GetHeight(), A.GetWidth()), a(&A)
1399{
1400 MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1401 "InverseMatrixCoefficient: "
1402 "Argument must be a square matrix.");
1403}
1404
1406{
1407 if (a) { a->SetTime(t); }
1409}
1410
1413 const IntegrationPoint &ip)
1414{
1415 a->Eval(M, T, ip);
1416 M.Invert();
1417}
1418
1420 : MatrixCoefficient(A.GetHeight(), A.GetWidth()), a(&A)
1421{
1422 MFEM_ASSERT(A.GetHeight() == A.GetWidth() && A.GetHeight() == 2,
1423 "ExponentialMatrixCoefficient: "
1424 << "Argument must be a square 2x2 matrix."
1425 << " Height = " << A.GetHeight()
1426 << ", Width = " << A.GetWidth());
1427}
1428
1430{
1431 if (a) { a->SetTime(t); }
1433}
1434
1437 const IntegrationPoint &ip)
1438{
1439 a->Eval(M, T, ip);
1440 M.Exponential();
1441}
1442
1445 : MatrixCoefficient(A.GetVDim(), B.GetVDim()), a(&A), b(&B),
1446 va(A.GetVDim()), vb(B.GetVDim())
1447{}
1448
1450{
1451 if (a) { a->SetTime(t); }
1452 if (b) { b->SetTime(t); }
1454}
1455
1457 const IntegrationPoint &ip)
1458{
1459 a->Eval(va, T, ip);
1460 b->Eval(vb, T, ip);
1461 M.SetSize(va.Size(), vb.Size());
1462 for (int i=0; i<va.Size(); i++)
1463 {
1464 for (int j=0; j<vb.Size(); j++)
1465 {
1466 M(i, j) = va[i] * vb[j];
1467 }
1468 }
1469}
1470
1472 : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(A), a(NULL), k(&K),
1473 vk(K.GetVDim())
1474{}
1475
1478 : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(0.0), a(&A), k(&K),
1479 vk(K.GetVDim())
1480{}
1481
1483{
1484 if (a) { a->SetTime(t); }
1485 if (k) { k->SetTime(t); }
1487}
1488
1490 const IntegrationPoint &ip)
1491{
1492 k->Eval(vk, T, ip);
1493 M.SetSize(vk.Size(), vk.Size());
1494 M = 0.0;
1495 real_t k2 = vk*vk;
1496 for (int i=0; i<vk.Size(); i++)
1497 {
1498 M(i, i) = k2;
1499 for (int j=0; j<vk.Size(); j++)
1500 {
1501 M(i, j) -= vk[i] * vk[j];
1502 }
1503 }
1504 M *= ((a == NULL ) ? aConst : a->Eval(T, ip) );
1505}
1506
1508 const IntegrationRule *irs[])
1509{
1510 real_t norm = 0.0;
1512
1513 for (int i = 0; i < mesh.GetNE(); i++)
1514 {
1515 tr = mesh.GetElementTransformation(i);
1516 const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1517 for (int j = 0; j < ir.GetNPoints(); j++)
1518 {
1519 const IntegrationPoint &ip = ir.IntPoint(j);
1520 tr->SetIntPoint(&ip);
1521 real_t val = fabs(coeff.Eval(*tr, ip));
1522 if (p < infinity())
1524 norm += ip.weight * tr->Weight() * pow(val, p);
1525 }
1526 else
1527 {
1528 if (norm < val)
1529 {
1530 norm = val;
1531 }
1532 }
1533 }
1534 }
1535 return norm;
1536}
1537
1539 const IntegrationRule *irs[])
1540{
1541 real_t norm = 0.0;
1543 int vdim = coeff.GetVDim();
1544 Vector vval(vdim);
1545 real_t val;
1546
1547 for (int i = 0; i < mesh.GetNE(); i++)
1548 {
1549 tr = mesh.GetElementTransformation(i);
1550 const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1551 for (int j = 0; j < ir.GetNPoints(); j++)
1552 {
1553 const IntegrationPoint &ip = ir.IntPoint(j);
1554 tr->SetIntPoint(&ip);
1555 coeff.Eval(vval, *tr, ip);
1556 if (p < infinity())
1557 {
1558 for (int idim(0); idim < vdim; ++idim)
1559 {
1560 norm += ip.weight * tr->Weight() * pow(fabs( vval(idim) ), p);
1561 }
1562 }
1563 else
1564 {
1565 for (int idim(0); idim < vdim; ++idim)
1566 {
1567 val = fabs(vval(idim));
1568 if (norm < val)
1569 {
1570 norm = val;
1571 }
1572 }
1573 }
1574 }
1575 }
1576
1577 return norm;
1578}
1579
1581 const IntegrationRule *irs[])
1582{
1583 real_t norm = LpNormLoop(p, coeff, mesh, irs);
1584
1585 if (p < infinity())
1586 {
1587 // negative quadrature weights may cause norm to be negative
1588 if (norm < 0.0)
1589 {
1590 norm = -pow(-norm, 1.0/p);
1591 }
1592 else
1593 {
1594 norm = pow(norm, 1.0/p);
1595 }
1596 }
1597
1598 return norm;
1599}
1600
1602 const IntegrationRule *irs[])
1603{
1604 real_t norm = LpNormLoop(p, coeff, mesh, irs);
1605
1606 if (p < infinity())
1607 {
1608 // negative quadrature weights may cause norm to be negative
1609 if (norm < 0.0)
1610 {
1611 norm = -pow(-norm, 1.0/p);
1612 }
1613 else
1614 {
1615 norm = pow(norm, 1.0/p);
1616 }
1617 }
1618
1619 return norm;
1620}
1621
1622#ifdef MFEM_USE_MPI
1624 const IntegrationRule *irs[])
1625{
1626 real_t loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1627 real_t glob_norm = 0;
1628
1629 MPI_Comm comm = pmesh.GetComm();
1630
1631 if (p < infinity())
1632 {
1633 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
1634 comm);
1635
1636 // negative quadrature weights may cause norm to be negative
1637 if (glob_norm < 0.0)
1638 {
1639 glob_norm = -pow(-glob_norm, 1.0/p);
1640 }
1641 else
1642 {
1643 glob_norm = pow(glob_norm, 1.0/p);
1644 }
1645 }
1646 else
1647 {
1648 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
1649 comm);
1650 }
1651
1652 return glob_norm;
1653}
1654
1656 const IntegrationRule *irs[])
1657{
1658 real_t loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1659 real_t glob_norm = 0;
1660
1661 MPI_Comm comm = pmesh.GetComm();
1662
1663 if (p < infinity())
1664 {
1665 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
1666 comm);
1667
1668 // negative quadrature weights may cause norm to be negative
1669 if (glob_norm < 0.0)
1670 {
1671 glob_norm = -pow(-glob_norm, 1.0/p);
1672 }
1673 else
1674 {
1675 glob_norm = pow(glob_norm, 1.0/p);
1676 }
1677 }
1678 else
1679 {
1680 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
1681 comm);
1682 }
1683
1684 return glob_norm;
1685}
1686#endif
1687
1691
1693{
1694 MFEM_VERIFY(index_ >= 0, "Index must be >= 0");
1695 MFEM_VERIFY(index_ < QuadF.GetVDim(),
1696 "Index must be < QuadratureFunction length");
1697 index = index_;
1698
1699 MFEM_VERIFY(length_ > 0, "Length must be > 0");
1700 MFEM_VERIFY(length_ <= QuadF.GetVDim() - index,
1701 "Length must be <= (QuadratureFunction length - index)");
1702
1703 vdim = length_;
1704}
1705
1708 const IntegrationPoint &ip)
1709{
1710 QuadF.HostRead();
1711
1712 const int el_idx = QuadF.GetSpace()->GetEntityIndex(T);
1713 // Handle the case of "interior boundary elements" and FaceQuadratureSpace
1714 // with FaceType::Boundary.
1715 if (el_idx < 0) { V = 0.0; return; }
1716
1717 const int ip_idx = QuadF.GetSpace()->GetPermutedIndex(el_idx, ip.index);
1718
1719 if (index == 0 && vdim == QuadF.GetVDim())
1720 {
1721 QuadF.GetValues(el_idx, ip_idx, V);
1722 }
1723 else
1724 {
1725 Vector temp;
1726 QuadF.GetValues(el_idx, ip_idx, temp);
1727 V.SetSize(vdim);
1728 for (int i = 0; i < vdim; i++)
1729 {
1730 V(i) = temp(index + i);
1731 }
1732 }
1733
1734 return;
1735}
1736
1741
1743 const QuadratureFunction &qf) : QuadF(qf)
1744{
1745 MFEM_VERIFY(qf.GetVDim() == 1, "QuadratureFunction's vdim must be 1");
1746}
1747
1749 const IntegrationPoint &ip)
1750{
1751 QuadF.HostRead();
1752 Vector temp(1);
1753 const int el_idx = QuadF.GetSpace()->GetEntityIndex(T);
1754 // Handle the case of "interior boundary elements" and FaceQuadratureSpace
1755 // with FaceType::Boundary.
1756 if (el_idx < 0) { return 0.0; }
1757 const int ip_idx = QuadF.GetSpace()->GetPermutedIndex(el_idx, ip.index);
1758 QuadF.GetValues(el_idx, ip_idx, temp);
1759 return temp[0];
1760}
1761
1763{
1764 qf = QuadF;
1765}
1766
1767
1770 : Vector(), storage(storage_), vdim(0), qs(qs_), qf(NULL)
1771{
1772 UseDevice(true);
1773}
1774
1777 CoefficientStorage storage_)
1778 : CoefficientVector(qs_, storage_)
1779{
1780 if (coeff == NULL)
1781 {
1782 SetConstant(1.0);
1783 }
1784 else
1785 {
1786 Project(*coeff);
1787 }
1788}
1789
1792 CoefficientStorage storage_)
1793 : CoefficientVector(qs_, storage_)
1794{
1795 Project(coeff);
1796}
1797
1800 CoefficientStorage storage_)
1801 : CoefficientVector(qs_, storage_)
1802{
1803 Project(coeff);
1804}
1805
1808 CoefficientStorage storage_)
1809 : CoefficientVector(qs_, storage_)
1810{
1811 Project(coeff);
1812}
1813
1815{
1816 vdim = 1;
1817 if (auto *const_coeff = dynamic_cast<ConstantCoefficient*>(&coeff))
1818 {
1819 SetConstant(const_coeff->constant);
1820 }
1821 else if (auto *qf_coeff = dynamic_cast<QuadratureFunctionCoefficient*>(&coeff))
1822 {
1823 MakeRef(qf_coeff->GetQuadFunction());
1824 }
1825 else
1826 {
1827 if (qf == nullptr) { qf = new QuadratureFunction(qs); }
1828 qf->SetVDim(1);
1829 coeff.Project(*qf);
1830 Vector::MakeRef(*qf, 0, qf->Size());
1831 }
1832}
1833
1835{
1836 vdim = coeff.GetVDim();
1837 if (auto *const_coeff = dynamic_cast<VectorConstantCoefficient*>(&coeff))
1838 {
1839 SetConstant(const_coeff->GetVec());
1840 }
1841 else if (auto *qf_coeff =
1842 dynamic_cast<VectorQuadratureFunctionCoefficient*>(&coeff))
1843 {
1844 MakeRef(qf_coeff->GetQuadFunction());
1845 }
1846 else
1847 {
1848 if (qf == nullptr) { qf = new QuadratureFunction(qs, vdim); }
1849 qf->SetVDim(vdim);
1850 coeff.Project(*qf);
1851 Vector::MakeRef(*qf, 0, qf->Size());
1852 }
1853}
1854
1856{
1857 if (auto *const_coeff = dynamic_cast<MatrixConstantCoefficient*>(&coeff))
1858 {
1859 SetConstant(const_coeff->GetMatrix());
1860 }
1861 else if (auto *const_sym_coeff =
1862 dynamic_cast<SymmetricMatrixConstantCoefficient*>(&coeff))
1863 {
1864 SetConstant(const_sym_coeff->GetMatrix());
1865 }
1866 else
1867 {
1868 auto *sym_coeff = dynamic_cast<SymmetricMatrixCoefficient*>(&coeff);
1869 const bool sym = sym_coeff && (storage & CoefficientStorage::SYMMETRIC);
1870 const int height = coeff.GetHeight();
1871 const int width = coeff.GetWidth();
1872 vdim = sym ? height*(height + 1)/2 : width*height;
1873
1874 if (qf == nullptr) { qf = new QuadratureFunction(qs, vdim); }
1875 qf->SetVDim(vdim);
1876 if (sym) { sym_coeff->ProjectSymmetric(*qf); }
1877 else { coeff.Project(*qf, transpose); }
1878 Vector::MakeRef(*qf, 0, qf->Size());
1879 }
1880}
1881
1883{
1884 Project(coeff, true);
1885}
1886
1888{
1889 vdim = qf_.GetVDim();
1890 const QuadratureSpaceBase *qs2 = qf_.GetSpace();
1891 MFEM_CONTRACT_VAR(qs2); // qs2 used only for asserts
1892 MFEM_VERIFY(qs2 != NULL, "Invalid QuadratureSpace.")
1893 MFEM_VERIFY(qs2->GetMesh() == qs.GetMesh(), "Meshes differ.");
1894 MFEM_VERIFY(qs2->GetOrder() == qs.GetOrder(), "Orders differ.");
1895 Vector::MakeRef(const_cast<QuadratureFunction&>(qf_), 0, qf_.Size());
1896}
1897
1899{
1900 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1901 vdim = 1;
1902 SetSize(nq);
1903 Vector::operator=(constant);
1904}
1905
1907{
1908 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1909 vdim = constant.Size();
1910 SetSize(nq*vdim);
1911 for (int iq = 0; iq < nq; ++iq)
1912 {
1913 for (int vd = 0; vd<vdim; ++vd)
1914 {
1915 (*this)[vd + iq*vdim] = constant[vd];
1916 }
1917 }
1918}
1919
1921{
1922 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1923 const int width = constant.Width();
1924 const int height = constant.Height();
1925 vdim = width*height;
1926 SetSize(nq*vdim);
1927 for (int iq = 0; iq < nq; ++iq)
1928 {
1929 for (int j = 0; j < width; ++j)
1930 {
1931 for (int i = 0; i < height; ++i)
1932 {
1933 (*this)[i + j*height + iq*vdim] = constant(i, j);
1934 }
1935 }
1936 }
1937}
1938
1940{
1941 const int nq = (storage & CoefficientStorage::CONSTANTS) ? 1 : qs.GetSize();
1942 const int height = constant.Height();
1943 const bool sym = storage & CoefficientStorage::SYMMETRIC;
1944 vdim = sym ? height*(height + 1)/2 : height*height;
1945 SetSize(nq*vdim);
1946 for (int iq = 0; iq < nq; ++iq)
1947 {
1948 for (int vd = 0; vd < vdim; ++vd)
1949 {
1950 const real_t value = sym ? constant.GetData()[vd] : constant(vd % height,
1951 vd / height);
1952 (*this)[vd + iq*vdim] = value;
1953 }
1954 }
1955}
1956
1957int CoefficientVector::GetVDim() const { return vdim; }
1958
1963
1964}
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
int Size() const
Return the logical size of the array.
Definition array.hpp:147
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
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) override
Fill the QuadratureFunction qf with the constant value.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
CrossCrossCoefficient(real_t A, VectorCoefficient &K)
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector curl coefficient at ip.
CurlGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
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:125
void Transpose()
(*this) = (*this)^t
void SetRow(int r, const real_t *row)
void GetColumnReference(int c, Vector &col)
Definition densemat.hpp:319
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:108
void Invert()
Replaces the current matrix with its inverse.
Definition densemat.cpp:707
real_t Trace() const
Trace of a square matrix.
Definition densemat.cpp:429
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:609
real_t Det() const
Definition densemat.cpp:516
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.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the determinant coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the scalar divergence coefficient at ip.
DivergenceGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
const Mesh * mesh
The Mesh object containing the element.
Definition eltrans.hpp:97
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:175
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
ExponentialMatrixCoefficient(MatrixCoefficient &A)
Construct the matrix coefficient. Result is .
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
std::function< real_t(const Vector &)> Function
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
std::function< real_t(const Vector &, real_t)> TDFunction
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the gradient vector coefficient at ip.
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...
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void Project(QuadratureFunction &qf) override
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:446
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
FiniteElementSpace * FESpace()
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition gridfunc.cpp:473
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:711
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
IsoparametricTransformation Transf
Definition eltrans.hpp:733
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition eltrans.cpp:587
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
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetTime(real_t t) override
Set the time for internally stored coefficients.
A standard isoparametric element transformation.
Definition eltrans.hpp:629
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition eltrans.hpp:668
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition eltrans.cpp:417
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) override
Set the time for internally stored coefficients.
void Eval(int i, Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
MatrixArrayVectorCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Set(int i, VectorCoefficient *c, bool own=true)
Set the coefficient located at the i-th row of the matrix. By this will take ownership of the Coeffic...
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.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip) override
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
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 Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Mesh data type.
Definition mesh.hpp:64
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition mesh.cpp:7635
const CoarseFineTransformations & GetRefinementTransforms() const
Definition mesh.cpp:11407
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
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
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
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
void SetTime(real_t t) override
Set the time for internally stored coefficients.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void UpdateCoefficient(int attr, Coefficient &coef)
Replace a single Coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient in the element described by T at the point ip.
QuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
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) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
ScalarMatrixProductCoefficient(real_t A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
ScalarVectorProductCoefficient(real_t A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Base class for symmetric matrix coefficients that optionally depend on time and space.
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 ...
DenseSymmetricMatrix mat_aux
Internal matrix used when evaluating this coefficient as a DenseMatrix.
A matrix coefficient that is constant in space and time.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the trace coefficient at ip.
TraceCoefficient(MatrixCoefficient &A)
Construct with the matrix.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set().
void SetTime(real_t t) override
Set the time for internally stored coefficients.
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 Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
VectorCrossProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two coefficients. Result is A x B.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
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...
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf.
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
VectorQuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
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 ...
void SetComponent(int index_, int length_)
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Vector data type.
Definition vector.hpp:82
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:498
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:931
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:147
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:506
Vector & operator=(const real_t *v)
Copy Size() entries from v.
Definition vector.cpp:196
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:622
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:236
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:47
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:548
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:391
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)
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:90
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:92
DenseTensor point_matrices[Geometry::NumGeom]
Definition ncmesh.hpp:96
Helper struct to convert a C++ type to an MPI type.