MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
coefficient.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 // Implementation of Coefficient class
13 
14 #include "fem.hpp"
15 
16 #include <cmath>
17 #include <limits>
18 
19 namespace mfem
20 {
21 
22 using 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  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;
35  const CoarseFineTransformations &cf = fine_mesh.GetRefinementTransforms();
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  const IntegrationPoint & ip)
53 {
54  int att = T.Attribute;
55  return (constants(att-1));
56 }
57 
58 void PWCoefficient::InitMap(const Array<int> & attr,
59  const Array<Coefficient*> & coefs)
60 {
61  MFEM_VERIFY(attr.Size() == coefs.Size(),
62  "PWCoefficient: "
63  "Attribute and coefficient arrays have incompatible "
64  "dimensions.");
65 
66  for (int i=0; i<attr.Size(); i++)
67  {
68  if (coefs[i] != NULL)
69  {
70  UpdateCoefficient(attr[i], *coefs[i]);
71  }
72  }
73 }
74 
76 {
78 
79  std::map<int, Coefficient*>::iterator p = pieces.begin();
80  for (; p != pieces.end(); p++)
81  {
82  if (p->second != NULL)
83  {
84  p->second->SetTime(t);
85  }
86  }
87 }
88 
90  const IntegrationPoint &ip)
91 {
92  const int att = T.Attribute;
93  std::map<int, Coefficient*>::const_iterator p = pieces.find(att);
94  if (p != pieces.end())
95  {
96  if ( p->second != NULL)
97  {
98  return p->second->Eval(T, ip);
99  }
100  }
101  return 0.0;
102 }
103 
105  const IntegrationPoint & ip)
106 {
107  double x[3];
108  Vector transip(x, 3);
109 
110  T.Transform(ip, transip);
111 
112  if (Function)
113  {
114  return Function(transip);
115  }
116  else
117  {
118  return TDFunction(transip, GetTime());
119  }
120 }
121 
123  const IntegrationPoint &ip)
124 {
125  Mesh *gf_mesh = GridF->FESpace()->GetMesh();
126  if (T.mesh == gf_mesh)
127  {
128  return GridF->GetValue(T, ip, Component);
129  }
130  else
131  {
132  IntegrationPoint coarse_ip;
133  ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
134  return GridF->GetValue(*coarse_T, coarse_ip, Component);
135  }
136 }
137 
139 {
140  if (Q1) { Q1->SetTime(t); }
141  if (Q2) { Q2->SetTime(t); }
142  this->Coefficient::SetTime(t);
143 }
144 
146  const IntegrationPoint &ip)
147 {
148  if (Q2)
149  {
150  return (*Transform2)(Q1->Eval(T, ip, GetTime()),
151  Q2->Eval(T, ip, GetTime()));
152  }
153  else
154  {
155  return (*Transform1)(Q1->Eval(T, ip, GetTime()));
156  }
157 }
158 
160 {
161  if (weight) { weight->SetTime(t); }
162  this->Coefficient::SetTime(t);
163 }
164 
166 {
167  MFEM_VERIFY(vcenter.Size() <= 3,
168  "SetDeltaCenter::Maximum number of dim supported is 3")
169  for (int i = 0; i < vcenter.Size(); i++) { center[i] = vcenter[i]; }
170  sdim = vcenter.Size();
171 }
172 
174 {
175  vcenter.SetSize(sdim);
176  vcenter = center;
177 }
178 
180  const IntegrationPoint &ip)
181 {
182  double w = Scale();
183  return weight ? weight->Eval(T, ip, GetTime())*w : w;
184 }
185 
187 {
188  if (c) { c->SetTime(t); }
189  this->Coefficient::SetTime(t);
190 }
191 
193  const IntegrationRule &ir)
194 {
195  Vector Mi;
196  M.SetSize(vdim, ir.GetNPoints());
197  for (int i = 0; i < ir.GetNPoints(); i++)
198  {
199  M.GetColumnReference(i, Mi);
200  const IntegrationPoint &ip = ir.IntPoint(i);
201  T.SetIntPoint(&ip);
202  Eval(Mi, T, ip);
203  }
204 }
205 
206 void PWVectorCoefficient::InitMap(const Array<int> & attr,
207  const Array<VectorCoefficient*> & coefs)
208 {
209  MFEM_VERIFY(attr.Size() == coefs.Size(),
210  "PWVectorCoefficient: "
211  "Attribute and coefficient arrays have incompatible "
212  "dimensions.");
213 
214  for (int i=0; i<attr.Size(); i++)
215  {
216  if (coefs[i] != NULL)
217  {
218  UpdateCoefficient(attr[i], *coefs[i]);
219  }
220  }
221 }
222 
224 {
225  MFEM_VERIFY(coef.GetVDim() == vdim,
226  "PWVectorCoefficient::UpdateCoefficient: "
227  "VectorCoefficient has incompatible dimension.");
228  pieces[attr] = &coef;
229 }
230 
232 {
234 
235  std::map<int, VectorCoefficient*>::iterator p = pieces.begin();
236  for (; p != pieces.end(); p++)
237  {
238  if (p->second != NULL)
239  {
240  p->second->SetTime(t);
241  }
242  }
243 }
244 
246  const IntegrationPoint &ip)
247 {
248  const int att = T.Attribute;
249  std::map<int, VectorCoefficient*>::const_iterator p = pieces.find(att);
250  if (p != pieces.end())
251  {
252  if ( p->second != NULL)
253  {
254  p->second->Eval(V, T, ip);
255  return;
256  }
257  }
258 
259  V.SetSize(vdim);
260  V = 0.0;
261 }
262 
264  const IntegrationPoint &ip)
265 {
266  double x[3];
267  Vector transip(x, 3);
268 
269  T.Transform(ip, transip);
270 
271  V.SetSize(vdim);
272  if (Function)
273  {
274  Function(transip, V);
275  }
276  else
277  {
278  TDFunction(transip, GetTime(), V);
279  }
280  if (Q)
281  {
282  V *= Q->Eval(T, ip, GetTime());
283  }
284 }
285 
287  : VectorCoefficient(dim), Coeff(dim), ownCoeff(dim)
288 {
289  for (int i = 0; i < dim; i++)
290  {
291  Coeff[i] = NULL;
292  ownCoeff[i] = true;
293  }
294 }
295 
297 {
298  for (int i = 0; i < vdim; i++)
299  {
300  if (Coeff[i]) { Coeff[i]->SetTime(t); }
301  }
303 }
304 
305 void VectorArrayCoefficient::Set(int i, Coefficient *c, bool own)
306 {
307  if (ownCoeff[i]) { delete Coeff[i]; }
308  Coeff[i] = c;
309  ownCoeff[i] = own;
310 }
311 
313 {
314  for (int i = 0; i < vdim; i++)
315  {
316  if (ownCoeff[i]) { delete Coeff[i]; }
317  }
318 }
319 
321  const IntegrationPoint &ip)
322 {
323  V.SetSize(vdim);
324  for (int i = 0; i < vdim; i++)
325  {
326  V(i) = this->Eval(i, T, ip);
327  }
328 }
329 
331  const GridFunction *gf)
332  : VectorCoefficient ((gf) ? gf -> VectorDim() : 0)
333 {
334  GridFunc = gf;
335 }
336 
338 {
339  GridFunc = gf; vdim = (gf) ? gf -> VectorDim() : 0;
340 }
341 
343  const IntegrationPoint &ip)
344 {
345  Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
346  if (T.mesh == gf_mesh)
347  {
348  GridFunc->GetVectorValue(T, ip, V);
349  }
350  else
351  {
352  IntegrationPoint coarse_ip;
353  ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
354  GridFunc->GetVectorValue(*coarse_T, coarse_ip, V);
355  }
356 }
357 
360 {
361  if (T.mesh == GridFunc->FESpace()->GetMesh())
362  {
363  GridFunc->GetVectorValues(T, ir, M);
364  }
365  else
366  {
367  VectorCoefficient::Eval(M, T, ir);
368  }
369 }
370 
372  const GridFunction *gf)
373  : VectorCoefficient((gf) ?
374  gf -> FESpace() -> GetMesh() -> SpaceDimension() : 0)
375 {
376  GridFunc = gf;
377 }
378 
380 {
381  GridFunc = gf; vdim = (gf) ?
382  gf -> FESpace() -> GetMesh() -> SpaceDimension() : 0;
383 }
384 
386  const IntegrationPoint &ip)
387 {
388  Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
389  if (T.mesh == gf_mesh)
390  {
391  GridFunc->GetGradient(T, V);
392  }
393  else
394  {
395  IntegrationPoint coarse_ip;
396  ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
397  GridFunc->GetGradient(*coarse_T, V);
398  }
399 }
400 
403 {
404  if (T.mesh == GridFunc->FESpace()->GetMesh())
405  {
406  GridFunc->GetGradients(T, ir, M);
407  }
408  else
409  {
410  VectorCoefficient::Eval(M, T, ir);
411  }
412 }
413 
415  const GridFunction *gf)
416  : VectorCoefficient(0)
417 {
418  SetGridFunction(gf);
419 }
420 
422 {
423  GridFunc = gf; vdim = (gf) ? gf -> CurlDim() : 0;
424 }
425 
427  const IntegrationPoint &ip)
428 {
429  Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
430  if (T.mesh == gf_mesh)
431  {
432  GridFunc->GetCurl(T, V);
433  }
434  else
435  {
436  IntegrationPoint coarse_ip;
437  ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
438  GridFunc->GetCurl(*coarse_T, V);
439  }
440 }
441 
443  const GridFunction *gf) : Coefficient()
444 {
445  GridFunc = gf;
446 }
447 
449  const IntegrationPoint &ip)
450 {
451  Mesh *gf_mesh = GridFunc->FESpace()->GetMesh();
452  if (T.mesh == gf_mesh)
453  {
454  return GridFunc->GetDivergence(T);
455  }
456  else
457  {
458  IntegrationPoint coarse_ip;
459  ElementTransformation *coarse_T = RefinedToCoarse(*gf_mesh, T, ip, coarse_ip);
460  return GridFunc->GetDivergence(*coarse_T);
461  }
462 }
463 
465 {
466  d.SetTime(t);
468 }
469 
471 {
472  dir = d_;
473  (*this).vdim = dir.Size();
474 }
475 
478 {
479  V = dir;
480  d.SetTime(GetTime());
481  V *= d.EvalDelta(T, ip);
482 }
483 
485 {
486  if (c) { c->SetTime(t); }
488 }
489 
491  const IntegrationPoint &ip)
492 {
493  V.SetSize(vdim);
494  if (active_attr[T.Attribute-1])
495  {
496  c->SetTime(GetTime());
497  c->Eval(V, T, ip);
498  }
499  else
500  {
501  V = 0.0;
502  }
503 }
504 
507 {
508  if (active_attr[T.Attribute-1])
509  {
510  c->SetTime(GetTime());
511  c->Eval(M, T, ir);
512  }
513  else
514  {
515  M.SetSize(vdim, ir.GetNPoints());
516  M = 0.0;
517  }
518 }
519 
520 void PWMatrixCoefficient::InitMap(const Array<int> & attr,
521  const Array<MatrixCoefficient*> & coefs)
522 {
523  MFEM_VERIFY(attr.Size() == coefs.Size(),
524  "PWMatrixCoefficient: "
525  "Attribute and coefficient arrays have incompatible "
526  "dimensions.");
527 
528  for (int i=0; i<attr.Size(); i++)
529  {
530  if (coefs[i] != NULL)
531  {
532  UpdateCoefficient(attr[i], *coefs[i]);
533  }
534  }
535 }
536 
538 {
539  MFEM_VERIFY(coef.GetHeight() == height,
540  "PWMatrixCoefficient::UpdateCoefficient: "
541  "MatrixCoefficient has incompatible height.");
542  MFEM_VERIFY(coef.GetWidth() == width,
543  "PWMatrixCoefficient::UpdateCoefficient: "
544  "MatrixCoefficient has incompatible width.");
545  if (symmetric)
546  {
547  MFEM_VERIFY(coef.IsSymmetric(),
548  "PWMatrixCoefficient::UpdateCoefficient: "
549  "MatrixCoefficient has incompatible symmetry.");
550  }
551  pieces[attr] = &coef;
552 }
553 
555 {
557 
558  std::map<int, MatrixCoefficient*>::iterator p = pieces.begin();
559  for (; p != pieces.end(); p++)
560  {
561  if (p->second != NULL)
562  {
563  p->second->SetTime(t);
564  }
565  }
566 }
567 
569  const IntegrationPoint &ip)
570 {
571  const int att = T.Attribute;
572  std::map<int, MatrixCoefficient*>::const_iterator p = pieces.find(att);
573  if (p != pieces.end())
574  {
575  if ( p->second != NULL)
576  {
577  p->second->Eval(K, T, ip);
578  return;
579  }
580  }
581 
582  K.SetSize(height, width);
583  K = 0.0;
584 }
585 
587 {
588  if (Q) { Q->SetTime(t); }
590 }
591 
593  const IntegrationPoint &ip)
594 {
595  double x[3];
596  Vector transip(x, 3);
597 
598  T.Transform(ip, transip);
599 
600  K.SetSize(height, width);
601 
602  if (symmetric) // Use SymmFunction (deprecated version)
603  {
604  MFEM_VERIFY(height == width && SymmFunction,
605  "MatrixFunctionCoefficient is not symmetric");
606 
607  Vector Ksym((width * (width + 1)) / 2); // 1x1: 1, 2x2: 3, 3x3: 6
608 
609  SymmFunction(transip, Ksym);
610 
611  // Copy upper triangular values from Ksym to the full matrix K
612  int os = 0;
613  for (int i=0; i<height; ++i)
614  {
615  for (int j=i; j<width; ++j)
616  {
617  const double Kij = Ksym[j - i + os];
618  K(i,j) = Kij;
619  if (j != i) { K(j,i) = Kij; }
620  }
621 
622  os += width - i;
623  }
624  }
625  else
626  {
627  if (Function)
628  {
629  Function(transip, K);
630  }
631  else if (TDFunction)
632  {
633  TDFunction(transip, GetTime(), K);
634  }
635  else
636  {
637  K = mat;
638  }
639  }
640 
641  if (Q)
642  {
643  K *= Q->Eval(T, ip, GetTime());
644  }
645 }
646 
649  const IntegrationPoint &ip)
650 {
651  MFEM_VERIFY(symmetric && height == width && SymmFunction,
652  "MatrixFunctionCoefficient is not symmetric");
653 
654  double x[3];
655  Vector transip(x, 3);
656 
657  T.Transform(ip, transip);
658 
659  K.SetSize((width * (width + 1)) / 2); // 1x1: 1, 2x2: 3, 3x3: 6
660 
661  if (SymmFunction)
662  {
663  SymmFunction(transip, K);
664  }
665 
666  if (Q)
667  {
668  K *= Q->Eval(T, ip, GetTime());
669  }
670 }
671 
673 {
674  if (Q) { Q->SetTime(t); }
676 }
677 
680  const IntegrationPoint &ip)
681 {
682  double x[3];
683  Vector transip(x, 3);
684 
685  T.Transform(ip, transip);
686 
687  K.SetSize(dim);
688 
689  if (Function)
690  {
691  Function(transip, K);
692  }
693  else if (TDFunction)
694  {
695  TDFunction(transip, GetTime(), K);
696  }
697  else
698  {
699  K = mat;
700  }
701 
702  if (Q)
703  {
704  K *= Q->Eval(T, ip, GetTime());
705  }
706 }
707 
709  : MatrixCoefficient (dim)
710 {
711  Coeff.SetSize(height*width);
712  ownCoeff.SetSize(height*width);
713  for (int i = 0; i < (height*width); i++)
714  {
715  Coeff[i] = NULL;
716  ownCoeff[i] = true;
717  }
718 }
719 
721 {
722  for (int i=0; i < height*width; i++)
723  {
724  if (Coeff[i]) { Coeff[i]->SetTime(t); }
725  }
727 }
728 
729 void MatrixArrayCoefficient::Set(int i, int j, Coefficient * c, bool own)
730 {
731  if (ownCoeff[i*width+j]) { delete Coeff[i*width+j]; }
732  Coeff[i*width+j] = c;
733  ownCoeff[i*width+j] = own;
734 }
735 
737 {
738  for (int i=0; i < height*width; i++)
739  {
740  if (ownCoeff[i]) { delete Coeff[i]; }
741  }
742 }
743 
745  const IntegrationPoint &ip)
746 {
747  K.SetSize(height, width);
748  for (int i = 0; i < height; i++)
749  {
750  for (int j = 0; j < width; j++)
751  {
752  K(i,j) = this->Eval(i, j, T, ip);
753  }
754  }
755 }
756 
758 {
759  if (c) { c->SetTime(t); }
761 }
762 
764  const IntegrationPoint &ip)
765 {
766  if (active_attr[T.Attribute-1])
767  {
768  c->SetTime(GetTime());
769  c->Eval(K, T, ip);
770  }
771  else
772  {
773  K.SetSize(height, width);
774  K = 0.0;
775  }
776 }
777 
779 {
780  if (a) { a->SetTime(t); }
781  if (b) { b->SetTime(t); }
782  this->Coefficient::SetTime(t);
783 }
784 
786 {
787  if (a) { a->SetTime(t); }
788  if (b) { b->SetTime(t); }
789  this->Coefficient::SetTime(t);
790 }
791 
793 {
794  if (a) { a->SetTime(t); }
795  if (b) { b->SetTime(t); }
796  this->Coefficient::SetTime(t);
797 }
798 
800 {
801  if (a) { a->SetTime(t); }
802  this->Coefficient::SetTime(t);
803 }
804 
807  : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
808 {
809  MFEM_ASSERT(A.GetVDim() == B.GetVDim(),
810  "InnerProductCoefficient: "
811  "Arguments have incompatible dimensions.");
812 }
813 
815 {
816  if (a) { a->SetTime(t); }
817  if (b) { b->SetTime(t); }
818  this->Coefficient::SetTime(t);
819 }
820 
822  const IntegrationPoint &ip)
823 {
824  a->Eval(va, T, ip);
825  b->Eval(vb, T, ip);
826  return va * vb;
827 }
828 
831  : a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
832 {
833  MFEM_ASSERT(A.GetVDim() == 2 && B.GetVDim() == 2,
834  "VectorRotProductCoefficient: "
835  "Arguments must have dimension equal to two.");
836 }
837 
839 {
840  if (a) { a->SetTime(t); }
841  if (b) { b->SetTime(t); }
842  this->Coefficient::SetTime(t);
843 }
844 
846  const IntegrationPoint &ip)
847 {
848  a->Eval(va, T, ip);
849  b->Eval(vb, T, ip);
850  return va[0] * vb[1] - va[1] * vb[0];
851 }
852 
854  : a(&A), ma(A.GetHeight(), A.GetWidth())
855 {
856  MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
857  "DeterminantCoefficient: "
858  "Argument must be a square matrix.");
859 }
860 
862 {
863  if (a) { a->SetTime(t); }
864  this->Coefficient::SetTime(t);
865 }
866 
868  const IntegrationPoint &ip)
869 {
870  a->Eval(ma, T, ip);
871  return ma.Det();
872 }
873 
875  : VectorCoefficient(dim),
876  ACoef(NULL), BCoef(NULL),
877  A(dim), B(dim),
878  alphaCoef(NULL), betaCoef(NULL),
879  alpha(1.0), beta(1.0)
880 {
881  A = 0.0; B = 0.0;
882 }
883 
885  VectorCoefficient &B_,
886  double alpha_, double beta_)
887  : VectorCoefficient(A_.GetVDim()),
888  ACoef(&A_), BCoef(&B_),
889  A(A_.GetVDim()), B(A_.GetVDim()),
890  alphaCoef(NULL), betaCoef(NULL),
891  alpha(alpha_), beta(beta_)
892 {
893  MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
894  "VectorSumCoefficient: "
895  "Arguments must have the same dimension.");
896 }
897 
899  VectorCoefficient &B_,
900  Coefficient &alpha_,
901  Coefficient &beta_)
902  : VectorCoefficient(A_.GetVDim()),
903  ACoef(&A_), BCoef(&B_),
904  A(A_.GetVDim()),
905  B(A_.GetVDim()),
906  alphaCoef(&alpha_),
907  betaCoef(&beta_),
908  alpha(0.0), beta(0.0)
909 {
910  MFEM_ASSERT(A_.GetVDim() == B_.GetVDim(),
911  "VectorSumCoefficient: "
912  "Arguments must have the same dimension.");
913 }
914 
916 {
917  if (ACoef) { ACoef->SetTime(t); }
918  if (BCoef) { BCoef->SetTime(t); }
919  if (alphaCoef) { alphaCoef->SetTime(t); }
920  if (betaCoef) { betaCoef->SetTime(t); }
922 }
923 
925  const IntegrationPoint &ip)
926 {
927  V.SetSize(A.Size());
928  if ( ACoef) { ACoef->Eval(A, T, ip); }
929  if ( BCoef) { BCoef->Eval(B, T, ip); }
930  if (alphaCoef) { alpha = alphaCoef->Eval(T, ip); }
931  if ( betaCoef) { beta = betaCoef->Eval(T, ip); }
932  add(alpha, A, beta, B, V);
933 }
934 
936  double A,
938  : VectorCoefficient(B.GetVDim()), aConst(A), a(NULL), b(&B)
939 {}
940 
942  Coefficient &A,
944  : VectorCoefficient(B.GetVDim()), aConst(0.0), a(&A), b(&B)
945 {}
946 
948 {
949  if (a) { a->SetTime(t); }
950  if (b) { b->SetTime(t); }
952 }
953 
955  const IntegrationPoint &ip)
956 {
957  double sa = (a == NULL) ? aConst : a->Eval(T, ip);
958  b->Eval(V, T, ip);
959  V *= sa;
960 }
961 
963  double tol_)
964  : VectorCoefficient(A.GetVDim()), a(&A), tol(tol_)
965 {}
966 
968 {
969  if (a) { a->SetTime(t); }
971 }
972 
974  const IntegrationPoint &ip)
975 {
976  a->Eval(V, T, ip);
977  double nv = V.Norml2();
978  V *= (nv > tol) ? (1.0/nv) : 0.0;
979 }
980 
984  : VectorCoefficient(3), a(&A), b(&B), va(A.GetVDim()), vb(B.GetVDim())
985 {
986  MFEM_ASSERT(A.GetVDim() == 3 && B.GetVDim() == 3,
987  "VectorCrossProductCoefficient: "
988  "Arguments must have dimension equal to three.");
989 }
990 
992 {
993  if (a) { a->SetTime(t); }
994  if (b) { b->SetTime(t); }
996 }
997 
999  const IntegrationPoint &ip)
1000 {
1001  a->Eval(va, T, ip);
1002  b->Eval(vb, T, ip);
1003  V.SetSize(3);
1004  V[0] = va[1] * vb[2] - va[2] * vb[1];
1005  V[1] = va[2] * vb[0] - va[0] * vb[2];
1006  V[2] = va[0] * vb[1] - va[1] * vb[0];
1007 }
1008 
1011  : VectorCoefficient(A.GetHeight()), a(&A), b(&B),
1012  ma(A.GetHeight(), A.GetWidth()), vb(B.GetVDim())
1013 {
1014  MFEM_ASSERT(A.GetWidth() == B.GetVDim(),
1015  "MatrixVectorProductCoefficient: "
1016  "Arguments have incompatible dimensions.");
1017 }
1018 
1020 {
1021  if (a) { a->SetTime(t); }
1022  if (b) { b->SetTime(t); }
1023  this->VectorCoefficient::SetTime(t);
1024 }
1025 
1027  const IntegrationPoint &ip)
1028 {
1029  a->Eval(ma, T, ip);
1030  b->Eval(vb, T, ip);
1031  V.SetSize(vdim);
1032  ma.Mult(vb, V);
1033 }
1034 
1036  const IntegrationPoint &ip)
1037 {
1038  M.SetSize(dim);
1039  M = 0.0;
1040  for (int d=0; d<dim; d++) { M(d,d) = 1.0; }
1041 }
1042 
1044  MatrixCoefficient &B,
1045  double alpha_, double beta_)
1046  : MatrixCoefficient(A.GetHeight(), A.GetWidth()),
1047  a(&A), b(&B), alpha(alpha_), beta(beta_),
1048  ma(A.GetHeight(), A.GetWidth())
1049 {
1050  MFEM_ASSERT(A.GetHeight() == B.GetHeight() && A.GetWidth() == B.GetWidth(),
1051  "MatrixSumCoefficient: "
1052  "Arguments must have the same dimensions.");
1053 }
1054 
1056 {
1057  if (a) { a->SetTime(t); }
1058  if (b) { b->SetTime(t); }
1059  this->MatrixCoefficient::SetTime(t);
1060 }
1061 
1063  const IntegrationPoint &ip)
1064 {
1065  b->Eval(M, T, ip);
1066  if ( beta != 1.0 ) { M *= beta; }
1067  a->Eval(ma, T, ip);
1068  M.Add(alpha, ma);
1069 }
1070 
1072  MatrixCoefficient &B)
1073  : MatrixCoefficient(A.GetHeight(), B.GetWidth()),
1074  a(&A), b(&B),
1075  ma(A.GetHeight(), A.GetWidth()),
1076  mb(B.GetHeight(), B.GetWidth())
1077 {
1078  MFEM_ASSERT(A.GetWidth() == B.GetHeight(),
1079  "MatrixProductCoefficient: "
1080  "Arguments must have compatible dimensions.");
1081 }
1082 
1084  const IntegrationPoint &ip)
1085 {
1086  a->Eval(ma, T, ip);
1087  b->Eval(mb, T, ip);
1088  Mult(ma, mb, M);
1089 }
1090 
1092  double A,
1093  MatrixCoefficient &B)
1094  : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(A), a(NULL), b(&B)
1095 {}
1096 
1098  Coefficient &A,
1099  MatrixCoefficient &B)
1100  : MatrixCoefficient(B.GetHeight(), B.GetWidth()), aConst(0.0), a(&A), b(&B)
1101 {}
1102 
1104 {
1105  if (a) { a->SetTime(t); }
1106  if (b) { b->SetTime(t); }
1107  this->MatrixCoefficient::SetTime(t);
1108 }
1109 
1112  const IntegrationPoint &ip)
1113 {
1114  double sa = (a == NULL) ? aConst : a->Eval(T, ip);
1115  b->Eval(M, T, ip);
1116  M *= sa;
1117 }
1118 
1120  : MatrixCoefficient(A.GetWidth(), A.GetHeight()), a(&A)
1121 {}
1122 
1124 {
1125  if (a) { a->SetTime(t); }
1126  this->MatrixCoefficient::SetTime(t);
1127 }
1128 
1131  const IntegrationPoint &ip)
1132 {
1133  a->Eval(M, T, ip);
1134  M.Transpose();
1135 }
1136 
1138  : MatrixCoefficient(A.GetHeight(), A.GetWidth()), a(&A)
1139 {
1140  MFEM_ASSERT(A.GetHeight() == A.GetWidth(),
1141  "InverseMatrixCoefficient: "
1142  "Argument must be a square matrix.");
1143 }
1144 
1146 {
1147  if (a) { a->SetTime(t); }
1148  this->MatrixCoefficient::SetTime(t);
1149 }
1150 
1153  const IntegrationPoint &ip)
1154 {
1155  a->Eval(M, T, ip);
1156  M.Invert();
1157 }
1158 
1160  VectorCoefficient &B)
1161  : MatrixCoefficient(A.GetVDim(), B.GetVDim()), a(&A), b(&B),
1162  va(A.GetVDim()), vb(B.GetVDim())
1163 {}
1164 
1166 {
1167  if (a) { a->SetTime(t); }
1168  if (b) { b->SetTime(t); }
1169  this->MatrixCoefficient::SetTime(t);
1170 }
1171 
1173  const IntegrationPoint &ip)
1174 {
1175  a->Eval(va, T, ip);
1176  b->Eval(vb, T, ip);
1177  M.SetSize(va.Size(), vb.Size());
1178  for (int i=0; i<va.Size(); i++)
1179  {
1180  for (int j=0; j<vb.Size(); j++)
1181  {
1182  M(i, j) = va[i] * vb[j];
1183  }
1184  }
1185 }
1186 
1188  : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(A), a(NULL), k(&K),
1189  vk(K.GetVDim())
1190 {}
1191 
1193  VectorCoefficient &K)
1194  : MatrixCoefficient(K.GetVDim(), K.GetVDim()), aConst(0.0), a(&A), k(&K),
1195  vk(K.GetVDim())
1196 {}
1197 
1199 {
1200  if (a) { a->SetTime(t); }
1201  if (k) { k->SetTime(t); }
1202  this->MatrixCoefficient::SetTime(t);
1203 }
1204 
1206  const IntegrationPoint &ip)
1207 {
1208  k->Eval(vk, T, ip);
1209  M.SetSize(vk.Size(), vk.Size());
1210  M = 0.0;
1211  double k2 = vk*vk;
1212  for (int i=0; i<vk.Size(); i++)
1213  {
1214  M(i, i) = k2;
1215  for (int j=0; j<vk.Size(); j++)
1216  {
1217  M(i, j) -= vk[i] * vk[j];
1218  }
1219  }
1220  M *= ((a == NULL ) ? aConst : a->Eval(T, ip) );
1221 }
1222 
1223 double LpNormLoop(double p, Coefficient &coeff, Mesh &mesh,
1224  const IntegrationRule *irs[])
1225 {
1226  double norm = 0.0;
1228 
1229  for (int i = 0; i < mesh.GetNE(); i++)
1230  {
1231  tr = mesh.GetElementTransformation(i);
1232  const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1233  for (int j = 0; j < ir.GetNPoints(); j++)
1234  {
1235  const IntegrationPoint &ip = ir.IntPoint(j);
1236  tr->SetIntPoint(&ip);
1237  double val = fabs(coeff.Eval(*tr, ip));
1238  if (p < infinity())
1239  {
1240  norm += ip.weight * tr->Weight() * pow(val, p);
1241  }
1242  else
1243  {
1244  if (norm < val)
1245  {
1246  norm = val;
1247  }
1248  }
1249  }
1250  }
1251  return norm;
1252 }
1253 
1254 double LpNormLoop(double p, VectorCoefficient &coeff, Mesh &mesh,
1255  const IntegrationRule *irs[])
1256 {
1257  double norm = 0.0;
1259  int vdim = coeff.GetVDim();
1260  Vector vval(vdim);
1261  double val;
1262 
1263  for (int i = 0; i < mesh.GetNE(); i++)
1264  {
1265  tr = mesh.GetElementTransformation(i);
1266  const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
1267  for (int j = 0; j < ir.GetNPoints(); j++)
1268  {
1269  const IntegrationPoint &ip = ir.IntPoint(j);
1270  tr->SetIntPoint(&ip);
1271  coeff.Eval(vval, *tr, ip);
1272  if (p < infinity())
1273  {
1274  for (int idim(0); idim < vdim; ++idim)
1275  {
1276  norm += ip.weight * tr->Weight() * pow(fabs( vval(idim) ), p);
1277  }
1278  }
1279  else
1280  {
1281  for (int idim(0); idim < vdim; ++idim)
1282  {
1283  val = fabs(vval(idim));
1284  if (norm < val)
1285  {
1286  norm = val;
1287  }
1288  }
1289  }
1290  }
1291  }
1292 
1293  return norm;
1294 }
1295 
1296 double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh,
1297  const IntegrationRule *irs[])
1298 {
1299  double norm = LpNormLoop(p, coeff, mesh, irs);
1300 
1301  if (p < infinity())
1302  {
1303  // negative quadrature weights may cause norm to be negative
1304  if (norm < 0.0)
1305  {
1306  norm = -pow(-norm, 1.0/p);
1307  }
1308  else
1309  {
1310  norm = pow(norm, 1.0/p);
1311  }
1312  }
1313 
1314  return norm;
1315 }
1316 
1317 double ComputeLpNorm(double p, VectorCoefficient &coeff, Mesh &mesh,
1318  const IntegrationRule *irs[])
1319 {
1320  double norm = LpNormLoop(p, coeff, mesh, irs);
1321 
1322  if (p < infinity())
1323  {
1324  // negative quadrature weights may cause norm to be negative
1325  if (norm < 0.0)
1326  {
1327  norm = -pow(-norm, 1.0/p);
1328  }
1329  else
1330  {
1331  norm = pow(norm, 1.0/p);
1332  }
1333  }
1334 
1335  return norm;
1336 }
1337 
1338 #ifdef MFEM_USE_MPI
1339 double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh,
1340  const IntegrationRule *irs[])
1341 {
1342  double loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1343  double glob_norm = 0;
1344 
1345  MPI_Comm comm = pmesh.GetComm();
1346 
1347  if (p < infinity())
1348  {
1349  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
1350 
1351  // negative quadrature weights may cause norm to be negative
1352  if (glob_norm < 0.0)
1353  {
1354  glob_norm = -pow(-glob_norm, 1.0/p);
1355  }
1356  else
1357  {
1358  glob_norm = pow(glob_norm, 1.0/p);
1359  }
1360  }
1361  else
1362  {
1363  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
1364  }
1365 
1366  return glob_norm;
1367 }
1368 
1369 double ComputeGlobalLpNorm(double p, VectorCoefficient &coeff, ParMesh &pmesh,
1370  const IntegrationRule *irs[])
1371 {
1372  double loc_norm = LpNormLoop(p, coeff, pmesh, irs);
1373  double glob_norm = 0;
1374 
1375  MPI_Comm comm = pmesh.GetComm();
1376 
1377  if (p < infinity())
1378  {
1379  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
1380 
1381  // negative quadrature weights may cause norm to be negative
1382  if (glob_norm < 0.0)
1383  {
1384  glob_norm = -pow(-glob_norm, 1.0/p);
1385  }
1386  else
1387  {
1388  glob_norm = pow(glob_norm, 1.0/p);
1389  }
1390  }
1391  else
1392  {
1393  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
1394  }
1395 
1396  return glob_norm;
1397 }
1398 #endif
1399 
1401  QuadratureFunction &qf)
1402  : VectorCoefficient(qf.GetVDim()), QuadF(qf), index(0) { }
1403 
1405 {
1406  MFEM_VERIFY(index_ >= 0, "Index must be >= 0");
1407  MFEM_VERIFY(index_ < QuadF.GetVDim(),
1408  "Index must be < QuadratureFunction length");
1409  index = index_;
1410 
1411  MFEM_VERIFY(length_ > 0, "Length must be > 0");
1412  MFEM_VERIFY(length_ <= QuadF.GetVDim() - index,
1413  "Length must be <= (QuadratureFunction length - index)");
1414 
1415  vdim = length_;
1416 }
1417 
1420  const IntegrationPoint &ip)
1421 {
1422  QuadF.HostRead();
1423 
1424  if (index == 0 && vdim == QuadF.GetVDim())
1425  {
1426  QuadF.GetElementValues(T.ElementNo, ip.index, V);
1427  }
1428  else
1429  {
1430  Vector temp;
1431  QuadF.GetElementValues(T.ElementNo, ip.index, temp);
1432  V.SetSize(vdim);
1433  for (int i = 0; i < vdim; i++)
1434  {
1435  V(i) = temp(index + i);
1436  }
1437  }
1438 
1439  return;
1440 }
1441 
1443  QuadratureFunction &qf) : QuadF(qf)
1444 {
1445  MFEM_VERIFY(qf.GetVDim() == 1, "QuadratureFunction's vdim must be 1");
1446 }
1447 
1449  const IntegrationPoint &ip)
1450 {
1451  QuadF.HostRead();
1452  Vector temp(1);
1453  QuadF.GetElementValues(T.ElementNo, ip.index, temp);
1454  return temp[0];
1455 }
1456 
1457 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
int GetVDim() const
Get the vector dimension.
Definition: gridfunc.hpp:825
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:703
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
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...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
Base class for vector Coefficients that optionally depend on time and space.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector curl coefficient at ip.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
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 ...
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.
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:465
NormalizedVectorCoefficient(VectorCoefficient &A, double tol=1e-6)
Return a vector normalized to a length of one.
void SetSize(int s)
Change the size of the DenseSymmetricMatrix to s x s.
Definition: symmat.cpp:32
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
double Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
double Det() const
Definition: densemat.cpp:436
void SetDeltaCenter(const Vector &center)
Set the center location of the delta function.
void SetTime(double t)
Set the time for internally stored coefficients.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
double GetTime()
Get the time for time dependent coefficients.
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:6197
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:792
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 ...
void SetTime(double t)
Set the time for internally stored coefficients.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
void SetTime(double t)
Set the time for internally stored coefficients.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector 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...
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
GradientGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a scalar grid function gf. The grid function is not owned by the coeff...
CrossCrossCoefficient(double A, VectorCoefficient &K)
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 coefficient at ip.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.hpp:48
void SetTime(double t)
Set the time for internally stored coefficients.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:300
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
ScalarMatrixProductCoefficient(double A, MatrixCoefficient &B)
Constructor with one coefficient. 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.
int GetHeight() const
Get the height of the matrix.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.cpp:75
double Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
void SetTime(double t)
Set the time for internally stored coefficients.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition: eltrans.hpp:401
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
DenseTensor point_matrices[Geometry::NumGeom]
Definition: ncmesh.hpp:77
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
double GetTime()
Get the time for time dependent coefficients.
ElementTransformation * RefinedToCoarse(Mesh &coarse_mesh, const ElementTransformation &T, const IntegrationPoint &ip, IntegrationPoint &coarse_ip)
Definition: coefficient.cpp:27
double LpNormLoop(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
QuadratureFunctionCoefficient(QuadratureFunction &qf)
Constructor with a quadrature function as input.
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:529
double b
Definition: lissajous.cpp:42
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:618
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:73
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
int GetVDim()
Returns dimension of the vector.
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
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 const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:442
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the determinant coefficient at ip.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetComponent(int index_, int length_)
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
bool IsSymmetric() const
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the scalar grid function.
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:1359
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition: eltrans.cpp:383
MPI_Comm GetComm() const
Definition: pmesh.hpp:288
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
VectorQuadratureFunctionCoefficient(QuadratureFunction &qf)
Constructor with a quadrature function as input.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Base class for Matrix Coefficients that optionally depend on time and space.
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:546
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
Definition: coefficient.cpp:51
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:285
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
double GetDivergence(ElementTransformation &tr) const
Definition: gridfunc.cpp:1456
virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9849
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
double a
Definition: lissajous.cpp:41
MatrixArrayCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the gradient vector coefficient at ip.
Class for integration point with weight.
Definition: intrules.hpp:25
void SetTime(double t)
Set the time for internally stored coefficients.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the scalar divergence coefficient at ip.
double GetTime()
Get the time for time dependent coefficients.
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
DivergenceGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(double t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
int dim
Definition: ex24.cpp:53
void SetDirection(const Vector &d_)
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: gridfunc.hpp:987
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf...
IsoparametricTransformation Transf
Definition: eltrans.hpp:464
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Definition: gridfunc.cpp:1709
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
void SetTime(double t)
Set the time for internally stored coefficients.
MatrixSumCoefficient(MatrixCoefficient &A, MatrixCoefficient &B, double alpha_=1.0, double beta_=1.0)
Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
RefCoord t[3]
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set()...
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
const double alpha
Definition: ex15.cpp:369
void GetCurl(ElementTransformation &tr, Vector &curl) const
Definition: gridfunc.cpp:1546
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:160
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
virtual void SetTime(double t)
Set the time for time dependent coefficients.
int GetWidth() const
Get the width of the matrix.
void SetTime(double t)
Set the time for internally stored coefficients.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:70
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
void SetTime(double t)
Set the time for internally stored coefficients.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
Definition: coefficient.cpp:89
void SetTime(double t)
Set the time for internally stored coefficients.
Class for parallel meshes.
Definition: pmesh.hpp:32
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetTime(double t)
Set the time for internally stored coefficients.
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:757
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
class Mesh * mesh
The Mesh object containing the element.
Definition: eltrans.hpp:84
void SetTime(double t)
Set the time for internally stored coefficients.
void GetGradient(ElementTransformation &tr, Vector &grad) const
Definition: gridfunc.cpp:1641
ScalarVectorProductCoefficient(double A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault...
void SetTime(double t)
Set the time for internally stored coefficients.