MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
coefficient.hpp
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#ifndef MFEM_COEFFICIENT
13#define MFEM_COEFFICIENT
14
15#include <functional>
16
17#include "../config/config.hpp"
18#include "../linalg/linalg.hpp"
19#include "intrules.hpp"
20#include "eltrans.hpp"
21
22namespace mfem
23{
24
25class Mesh;
26class QuadratureSpaceBase;
27class QuadratureFunction;
28
29#ifdef MFEM_USE_MPI
30class ParMesh;
31#endif
32
33
34/** @brief Base class Coefficients that optionally depend on space and time.
35 These are used by the BilinearFormIntegrator, LinearFormIntegrator, and
36 NonlinearFormIntegrator classes to represent the physical coefficients in
37 the PDEs that are being discretized. This class can also be used in a more
38 general way to represent functions that don't necessarily belong to a FE
39 space, e.g., to project onto GridFunctions to use as initial conditions,
40 exact solutions, etc. See, e.g., ex4 or ex22 for these uses. */
42{
43protected:
45
46public:
47 Coefficient() { time = 0.; }
48
49 /// Set the time for time dependent coefficients
50 virtual void SetTime(real_t t) { time = t; }
51
52 /// Get the time for time dependent coefficients
53 real_t GetTime() { return time; }
54
55 /** @brief Evaluate the coefficient in the element described by @a T at the
56 point @a ip. */
57 /** @note When this method is called, the caller must make sure that the
58 IntegrationPoint associated with @a T is the same as @a ip. This can be
59 achieved by calling T.SetIntPoint(&ip). */
61 const IntegrationPoint &ip) = 0;
62
63 /** @brief Evaluate the coefficient in the element described by @a T at the
64 point @a ip at time @a t. */
65 /** @note When this method is called, the caller must make sure that the
66 IntegrationPoint associated with @a T is the same as @a ip. This can be
67 achieved by calling T.SetIntPoint(&ip). */
69 const IntegrationPoint &ip, real_t t)
70 {
71 SetTime(t);
72 return Eval(T, ip);
73 }
74
75 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
76 /// the quadrature points.
77 virtual void Project(QuadratureFunction &qf);
78
79 virtual ~Coefficient() { }
80};
81
82
83/// A coefficient that is constant across space and time
85{
86public:
88
89 /// c is value of constant function
90 explicit ConstantCoefficient(real_t c = 1.0) { constant=c; }
91
92 /// Evaluate the coefficient at @a ip.
94 const IntegrationPoint &ip) override
95 { return (constant); }
96
97 /// Fill the QuadratureFunction @a qf with the constant value.
98 void Project(QuadratureFunction &qf) override;
99};
100
101/** @brief A piecewise constant coefficient with the constants keyed
102 off the element attribute numbers. */
104{
105private:
106 Vector constants;
107
108public:
109
110 /// Constructs a piecewise constant coefficient in NumOfSubD subdomains
111 explicit PWConstCoefficient(int NumOfSubD = 0) : constants(NumOfSubD)
112 { constants = 0.0; }
113
114 /// Construct the constant coefficient using a vector of constants.
115 /** @a c should be a vector defined by attributes, so for region with
116 attribute @a i @a c[i-1] is the coefficient in that region */
118 { constants.SetSize(c.Size()); constants=c; }
119
120 /// Update the constants with vector @a c.
121 void UpdateConstants(Vector &c) { constants.SetSize(c.Size()); constants=c; }
122
123 /// Return a reference to the i-th constant
124 real_t &operator()(int i) { return constants(i-1); }
125
126 /// Set the constants for all attributes to constant @a c.
127 void operator=(real_t c) { constants = c; }
128
129 /// Returns the number of constants representing different attributes.
130 int GetNConst() { return constants.Size(); }
131
132 /// Evaluate the coefficient.
134 const IntegrationPoint &ip) override;
135};
136
137/** @brief A piecewise coefficient with the pieces keyed off the element
138 attribute numbers.
139
140 A value of zero will be returned for any missing attribute numbers.
141
142 This object will not assume ownership of any Coefficient objects
143 passed to it. Consequently, the caller must ensure that the
144 individual Coefficient objects are not deleted while this
145 PWCoefficient is still in use.
146
147 \note The keys may either be domain attribute numbers or boundary
148 attribute numbers. If the PWCoefficient is used with a domain
149 integrator the keys are assumed to be domain attribute
150 numbers. Similarly, if the PWCoefficient is used with a boundary
151 integrator the keys are assumed to be boundary attribute numbers.
152*/
154{
155private:
156 /** Internal data structure to store pointers to the appropriate
157 coefficients for different regions of the mesh. The keys used
158 in the map are the mesh attribute numbers (either element
159 attribute or boundary element attribute depending upon
160 context). The values returned for any missing attributes will
161 be zero. The coefficient pointers may be NULL in which case a
162 value of zero is returned.
163
164 The Coefficient objects contained in this map are NOT owned by
165 this PWCoefficient object. This means that they will not be
166 deleted when this object is deleted also the caller must ensure
167 that the various Coefficient objects are not deleted while this
168 PWCoefficient is still needed.
169 */
170 std::map<int, Coefficient*> pieces;
171
172 /** Convenience function to check for compatible array lengths,
173 loop over the arrays, and add their attribute/Coefficient pairs
174 to the internal data structure.
175 */
176 void InitMap(const Array<int> & attr,
177 const Array<Coefficient*> & coefs);
178
179public:
180
181 /// Constructs a piecewise coefficient
182 explicit PWCoefficient() {}
183
184 /// Construct the coefficient using arrays describing the pieces
185 /** \param attr - an array of attribute numbers for each piece
186 \param coefs - the corresponding array of Coefficient pointers
187 Any missing attributes or NULL coefficient pointers will result in a
188 value of zero being returned for that attribute.
189
190 \note Ownership of the Coefficient objects will NOT be
191 transferred to this object.
192 */
194 const Array<Coefficient*> & coefs)
195 { InitMap(attr, coefs); }
196
197 /// Set the time for time dependent coefficients
198 void SetTime(real_t t) override;
199
200 /// Replace a set of coefficients
202 const Array<Coefficient*> & coefs)
203 { InitMap(attr, coefs); }
204
205 /// Replace a single Coefficient for a particular attribute
206 void UpdateCoefficient(int attr, Coefficient & coef)
207 { pieces[attr] = &coef; }
208
209 /// Remove a single Coefficient for a particular attribute
210 void ZeroCoefficient(int attr)
211 { pieces.erase(attr); }
212
213 /// Evaluate the coefficient.
215 const IntegrationPoint &ip) override;
216};
217
218/// A general function coefficient
220{
221protected:
222 std::function<real_t(const Vector &)> Function;
223 std::function<real_t(const Vector &, real_t)> TDFunction;
224
225public:
226 /// Define a time-independent coefficient from a std function
227 /** \param F time-independent std::function */
228 FunctionCoefficient(std::function<real_t(const Vector &)> F)
229 : Function(std::move(F))
230 { }
231
232 /// Define a time-dependent coefficient from a std function
233 /** \param TDF time-dependent function */
234 FunctionCoefficient(std::function<real_t(const Vector &, real_t)> TDF)
235 : TDFunction(std::move(TDF))
236 { }
237
238 /// (DEPRECATED) Define a time-independent coefficient from a C-function
239 /** @deprecated Use the method where the C-function, @a f, uses a const
240 Vector argument instead of Vector. */
241 MFEM_DEPRECATED FunctionCoefficient(real_t (*f)(Vector &))
242 {
243 Function = reinterpret_cast<real_t(*)(const Vector&)>(f);
244 TDFunction = NULL;
245 }
246
247 /// (DEPRECATED) Define a time-dependent coefficient from a C-function
248 /** @deprecated Use the method where the C-function, @a tdf, uses a const
249 Vector argument instead of Vector. */
250 MFEM_DEPRECATED FunctionCoefficient(real_t (*tdf)(Vector &, real_t))
251 {
252 Function = NULL;
253 TDFunction = reinterpret_cast<real_t(*)(const Vector&,real_t)>(tdf);
254 }
255
256 /// Evaluate the coefficient at @a ip.
258 const IntegrationPoint &ip) override;
259};
260
261/// A common base class for returning individual components of the domain's
262/// Cartesian coordinates.
264{
265protected:
266 int comp;
268
269 /// @a comp_ index of the desired component (0 -> x, 1 -> y, 2 -> z)
270 CartesianCoefficient(int comp_) : comp(comp_), transip(3) {}
271
272public:
273 /// Evaluate the coefficient at @a ip.
275 const IntegrationPoint &ip) override;
276};
277
278/// Scalar coefficient which returns the x-component of the evaluation point
284
285/// Scalar coefficient which returns the y-component of the evaluation point
291
292/// Scalar coefficient which returns the z-component of the evaluation point
298
299/// Scalar coefficient which returns the radial distance from the axis of
300/// the evaluation point in the cylindrical coordinate system
302{
303private:
304 mutable Vector transip;
305
306public:
308
309 /// Evaluate the coefficient at @a ip.
311 const IntegrationPoint &ip) override;
312};
313
314/// Scalar coefficient which returns the angular position or azimuth (often
315/// denoted by theta) of the evaluation point in the cylindrical coordinate
316/// system
318{
319private:
320 mutable Vector transip;
321
322public:
324
325 /// Evaluate the coefficient at @a ip.
327 const IntegrationPoint &ip) override;
328};
329
330/// Scalar coefficient which returns the height or altitude of
331/// the evaluation point in the cylindrical coordinate system
333
334/// Scalar coefficient which returns the radial distance from the origin of
335/// the evaluation point in the spherical coordinate system
337{
338private:
339 mutable Vector transip;
340
341public:
343
344 /// Evaluate the coefficient at @a ip.
346 const IntegrationPoint &ip) override;
347};
348
349/// Scalar coefficient which returns the azimuthal angle (often denoted by phi)
350/// of the evaluation point in the spherical coordinate system
352{
353private:
354 mutable Vector transip;
355
356public:
358
359 /// Evaluate the coefficient at @a ip.
361 const IntegrationPoint &ip) override;
362};
363
364/// Scalar coefficient which returns the polar angle (often denoted by theta)
365/// of the evaluation point in the spherical coordinate system
367{
368private:
369 mutable Vector transip;
370
371public:
373
374 /// Evaluate the coefficient at @a ip.
376 const IntegrationPoint &ip) override;
377};
378
379class GridFunction;
380
381/// Coefficient defined by a GridFunction. This coefficient is mesh dependent.
383{
384private:
385 const GridFunction *GridF;
386 int Component;
387
388public:
389 GridFunctionCoefficient() : GridF(NULL), Component(1) { }
390 /** Construct GridFunctionCoefficient from a given GridFunction, and
391 optionally specify a component to use if it is a vector GridFunction. */
392 GridFunctionCoefficient (const GridFunction *gf, int comp = 1)
393 { GridF = gf; Component = comp; }
394
395 /// Set the internal GridFunction
396 void SetGridFunction(const GridFunction *gf) { GridF = gf; }
397
398 /// Get the internal GridFunction
399 const GridFunction * GetGridFunction() const { return GridF; }
400
401 /// Evaluate the coefficient at @a ip.
403 const IntegrationPoint &ip) override;
404
405 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
406 /// the quadrature points.
407 ///
408 /// This function uses the efficient QuadratureFunction::ProjectGridFunction
409 /// to fill the QuadratureFunction.
410 void Project(QuadratureFunction &qf) override;
411};
412
413
414/** @brief A coefficient that depends on 1 or 2 parent coefficients and a
415 transformation rule represented by a C-function.
416
417 $ C(x,t) = T(Q1(x,t)) $ or $ C(x,t) = T(Q1(x,t), Q2(x,t)) $
418
419 where T is the transformation rule, and Q1/Q2 are the parent coefficients.*/
421{
422private:
423 Coefficient * Q1;
424 Coefficient * Q2;
425 std::function<real_t(real_t)> Transform1;
426 std::function<real_t(real_t, real_t)> Transform2;
427
428public:
430 : Q1(q), Transform1(std::move(F)) { Q2 = 0; Transform2 = 0; }
432 std::function<real_t(real_t, real_t)> F)
433 : Q1(q1), Q2(q2), Transform2(std::move(F)) { Transform1 = 0; }
434
435 /// Set the time for internally stored coefficients
436 void SetTime(real_t t) override;
437
438 /// Evaluate the coefficient at @a ip.
439 real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override;
440};
441
442/** @brief Delta function coefficient optionally multiplied by a weight
443 coefficient and a scaled time dependent C-function.
444
445 $ F(x,t) = w(x,t) s T(t) d(x - xc) $
446
447 where w is the optional weight coefficient, @a s is a scale factor
448 T is an optional time-dependent function and d is a delta function.
449
450 WARNING this cannot be used as a normal coefficient. The usual Eval
451 method is disabled. */
453{
454protected:
457 int sdim;
459
460public:
461
462 /// Construct a unit delta function centered at (0.0,0.0,0.0)
464 {
465 center[0] = center[1] = center[2] = 0.; scale = 1.; tol = 1e-12;
466 weight = NULL; sdim = 0; tdf = NULL;
467 }
468
469 /// Construct a delta function scaled by @a s and centered at (x,0.0,0.0)
471 {
472 center[0] = x; center[1] = 0.; center[2] = 0.; scale = s; tol = 1e-12;
473 weight = NULL; sdim = 1; tdf = NULL;
474 }
475
476 /// Construct a delta function scaled by @a s and centered at (x,y,0.0)
478 {
479 center[0] = x; center[1] = y; center[2] = 0.; scale = s; tol = 1e-12;
480 weight = NULL; sdim = 2; tdf = NULL;
481 }
482
483 /// Construct a delta function scaled by @a s and centered at (x,y,z)
485 {
486 center[0] = x; center[1] = y; center[2] = z; scale = s; tol = 1e-12;
487 weight = NULL; sdim = 3; tdf = NULL;
488 }
489
490 /// Set the time for internally stored coefficients
491 void SetTime(real_t t) override;
492
493 /// Set the center location of the delta function.
494 void SetDeltaCenter(const Vector& center);
495
496 /// Set the scale value multiplying the delta function.
497 void SetScale(real_t s_) { scale = s_; }
498
499 /// Set a time-dependent function that multiplies the Scale().
500 void SetFunction(real_t (*f)(real_t)) { tdf = f; }
501
502 /** @brief Set the tolerance used during projection onto GridFunction to
503 identify the Mesh vertex where the Center() of the delta function
504 lies. (default 1e-12)*/
505 void SetTol(real_t tol_) { tol = tol_; }
506
507 /// Set a weight Coefficient that multiplies the DeltaCoefficient.
508 /** The weight Coefficient multiplies the value returned by EvalDelta() but
509 not the value returned by Scale().
510 The weight Coefficient is also used as the L2-weight function when
511 projecting the DeltaCoefficient onto a GridFunction, so that the weighted
512 integral of the projection is exactly equal to the Scale(). */
513 void SetWeight(Coefficient *w) { weight = w; }
514
515 /// Return a pointer to a c-array representing the center of the delta
516 /// function.
517 const real_t *Center() { return center; }
518
519 /** @brief Return the scale factor times the optional time dependent
520 function. Returns $ s T(t) $ with $ T(t) = 1 $ when
521 not set by the user. */
522 real_t Scale() { return tdf ? (*tdf)(GetTime())*scale : scale; }
523
524 /// Return the tolerance used to identify the mesh vertices
525 real_t Tol() { return tol; }
526
527 /// See SetWeight() for description of the weight Coefficient.
528 Coefficient *Weight() { return weight; }
529
530 /// Write the center of the delta function into @a center.
532
533 /// The value of the function assuming we are evaluating at the delta center.
535 /** @brief A DeltaFunction cannot be evaluated. Calling this method will
536 cause an MFEM error, terminating the application. */
538 { mfem_error("DeltaCoefficient::Eval"); return 0.; }
539 virtual ~DeltaCoefficient() { delete weight; }
540};
541
542/** @brief Derived coefficient that takes the value of the parent coefficient
543 for the active attributes and is zero otherwise. */
545{
546private:
547 Coefficient *c;
548 Array<int> active_attr;
549
550public:
551 /** @brief Construct with a parent coefficient and an array with
552 ones marking the attributes on which this coefficient should be
553 active. */
555 { c = &c_; attr.Copy(active_attr); }
556
557 /// Set the time for internally stored coefficients
558 void SetTime(real_t t) override;
559
560 /// Evaluate the coefficient at @a ip.
562 { return active_attr[T.Attribute-1] ? c->Eval(T, ip, GetTime()) : 0.0; }
563};
564
565/// Base class for vector Coefficients that optionally depend on time and space.
567{
568protected:
569 int vdim;
571
572public:
573 /// Initialize the VectorCoefficient with vector dimension @a vd.
574 VectorCoefficient(int vd) { vdim = vd; time = 0.; }
575
576 /// Set the time for time dependent coefficients
577 virtual void SetTime(real_t t) { time = t; }
578
579 /// Get the time for time dependent coefficients
580 real_t GetTime() { return time; }
581
582 /// Returns dimension of the vector.
583 int GetVDim() { return vdim; }
584
585 /** @brief Evaluate the vector coefficient in the element described by @a T
586 at the point @a ip, storing the result in @a V. */
587 /** @note When this method is called, the caller must make sure that the
588 IntegrationPoint associated with @a T is the same as @a ip. This can be
589 achieved by calling T.SetIntPoint(&ip). */
590 virtual void Eval(Vector &V, ElementTransformation &T,
591 const IntegrationPoint &ip) = 0;
592
593 /** @brief Evaluate the vector coefficient in the element described by @a T
594 at all points of @a ir, storing the result in @a M. */
595 /** The dimensions of @a M are GetVDim() by ir.GetNPoints() and they must be
596 set by the implementation of this method.
597
598 The general implementation provided by the base class (using the Eval
599 method for one IntegrationPoint at a time) can be overloaded for more
600 efficient implementation.
601
602 @note The IntegrationPoint associated with @a T is not used, and this
603 method will generally modify this IntegrationPoint associated with @a T.
604 */
605 virtual void Eval(DenseMatrix &M, ElementTransformation &T,
606 const IntegrationRule &ir);
607
608 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
609 /// the quadrature points.
610 ///
611 /// The @a vdim of the VectorCoefficient should be equal to the @a vdim of
612 /// the QuadratureFunction.
613 virtual void Project(QuadratureFunction &qf);
614
615 virtual ~VectorCoefficient() { }
616};
617
618
619/// Vector coefficient that is constant in space and time.
621{
622private:
623 Vector vec;
624public:
625 /// Construct the coefficient with constant vector @a v.
627 : VectorCoefficient(v.Size()), vec(v) { }
629
630 /// Evaluate the vector coefficient at @a ip.
632 const IntegrationPoint &ip) override { V = vec; }
633
634 /// Return a reference to the constant vector in this class.
635 const Vector& GetVec() const { return vec; }
636};
637
638/** @brief A piecewise vector-valued coefficient with the pieces keyed off the
639 element attribute numbers.
640
641 A value of zero will be returned for any missing attribute numbers.
642
643 This object will not assume ownership of any VectorCoefficient
644 objects passed to it. Consequently, the caller must ensure that
645 the individual VectorCoefficient objects are not deleted while
646 this PWVectorCoefficient is still in use.
647
648 \note The keys may either be domain attribute numbers or boundary
649 attribute numbers. If the PWVectorCoefficient is used with a
650 domain integrator the keys are assumed to be domain attribute
651 numbers. Similarly, if the PWVectorCoefficient is used with a
652 boundary integrator the keys are assumed to be boundary attribute
653 numbers.
654*/
656{
657private:
658 /** Internal data structure to store pointers to the appropriate
659 coefficients for different regions of the mesh. The keys used
660 in the map are the mesh attribute numbers (either element
661 attribute or boundary element attribute depending upon
662 context). The values returned for any missing attributes will
663 be zero. The coefficient pointers may be NULL in which case a
664 value of zero is returned.
665
666 The VectorCoefficient objects contained in this map are NOT
667 owned by this PWVectorCoefficient object. This means that they
668 will not be deleted when this object is deleted also the caller
669 must ensure that the various VectorCoefficient objects are not
670 deleted while this PWVectorCoefficient is still needed.
671 */
672 std::map<int, VectorCoefficient*> pieces;
673
674 /** Convenience function to check for compatible array lengths,
675 loop over the arrays, and add their attribute/VectorCoefficient
676 pairs to the internal data structure.
677 */
678 void InitMap(const Array<int> & attr,
679 const Array<VectorCoefficient*> & coefs);
680
681public:
682
683 /// Constructs a piecewise vector coefficient of dimension vd
684 explicit PWVectorCoefficient(int vd): VectorCoefficient(vd) {}
685
686 /// Construct the coefficient using arrays describing the pieces
687 /** \param vd - dimension of the vector-valued result
688 \param attr - an array of attribute numbers for each piece
689 \param coefs - the corresponding array of VectorCoefficient pointers
690 Any missing attributes or NULL coefficient pointers will result in a
691 zero vector being returned for that attribute.
692
693 \note Ownership of the VectorCoefficient objects will NOT be
694 transferred to this object.
695 */
696 PWVectorCoefficient(int vd, const Array<int> & attr,
697 const Array<VectorCoefficient*> & coefs)
698 : VectorCoefficient(vd) { InitMap(attr, coefs); }
699
700 /// Set the time for time dependent coefficients
701 void SetTime(real_t t) override;
702
703 /// Replace a set of coefficients
705 const Array<VectorCoefficient*> & coefs)
706 { InitMap(attr, coefs); }
707
708 /// Replace a single Coefficient for a particular attribute
709 void UpdateCoefficient(int attr, VectorCoefficient & coef);
710
711 /// Remove a single VectorCoefficient for a particular attribute
712 void ZeroCoefficient(int attr)
713 { pieces.erase(attr); }
714
715 /// Evaluate the coefficient.
717 const IntegrationPoint &ip) override;
719};
720
721/// A vector coefficient which returns the physical location of the
722/// evaluation point in the Cartesian coordinate system.
724{
725public:
726
728
730 /// Evaluate the vector coefficient at @a ip.
732 const IntegrationPoint &ip) override;
733
735};
736
737/// A general vector function coefficient
739{
740private:
741 std::function<void(const Vector &, Vector &)> Function;
742 std::function<void(const Vector &, real_t, Vector &)> TDFunction;
743 Coefficient *Q;
744
745public:
746 /// Define a time-independent vector coefficient from a std function
747 /** \param dim - the size of the vector
748 \param F - time-independent function
749 \param q - optional scalar Coefficient to scale the vector coefficient */
751 std::function<void(const Vector &, Vector &)> F,
752 Coefficient *q = nullptr)
753 : VectorCoefficient(dim), Function(std::move(F)), Q(q)
754 { }
755
756 /// Define a time-dependent vector coefficient from a std function
757 /** \param dim - the size of the vector
758 \param TDF - time-dependent function
759 \param q - optional scalar Coefficient to scale the vector coefficient */
761 std::function<void(const Vector &, real_t, Vector &)> TDF,
762 Coefficient *q = nullptr)
763 : VectorCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
764 { }
765
767 /// Evaluate the vector coefficient at @a ip.
769 const IntegrationPoint &ip) override;
770
772};
773
774/** @brief Vector coefficient defined by an array of scalar coefficients.
775 Coefficients that are not set will evaluate to zero in the vector. This
776 object takes ownership of the array of coefficients inside it and deletes
777 them at object destruction. */
779{
780private:
782 Array<bool> ownCoeff;
783
784public:
785 /** @brief Construct vector of dim coefficients. The actual coefficients
786 still need to be added with Set(). */
787 explicit VectorArrayCoefficient(int dim);
788
789 /// Set the time for internally stored coefficients
790 void SetTime(real_t t) override;
791
792 /// Returns i'th coefficient.
793 Coefficient* GetCoeff(int i) { return Coeff[i]; }
794
795 /// Returns the entire array of coefficients.
796 Coefficient **GetCoeffs() { return Coeff; }
797
798 /// Sets coefficient in the vector.
799 void Set(int i, Coefficient *c, bool own=true);
800
801 /// Set ownership of the i'th coefficient
802 void SetOwnership(int i, bool own) { ownCoeff[i] = own; }
803
804 /// Get ownership of the i'th coefficient
805 bool GetOwnership(int i) const { return ownCoeff[i]; }
806
807 /// Evaluates i'th component of the vector of coefficients and returns the
808 /// value.
810 { return Coeff[i] ? Coeff[i]->Eval(T, ip, GetTime()) : 0.0; }
811
813 /** @brief Evaluate the coefficient. Each element of vector V comes from the
814 associated array of scalar coefficients. */
816 const IntegrationPoint &ip) override;
817
818 /// Destroys vector coefficient.
819 virtual ~VectorArrayCoefficient();
820};
821
822/// Vector coefficient defined by a vector GridFunction
824{
825protected:
827
828public:
829 /** @brief Construct an empty coefficient. Calling Eval() before the grid
830 function is set will cause a segfault. */
832
833 /** @brief Construct the coefficient with grid function @a gf. The
834 grid function is not owned by the coefficient. */
836
837 /** @brief Set the grid function for this coefficient. Also sets the Vector
838 dimension to match that of the @a gf. */
839 void SetGridFunction(const GridFunction *gf);
840
841 /// Returns a pointer to the grid function in this Coefficient
842 const GridFunction * GetGridFunction() const { return GridFunc; }
843
844 /// Evaluate the vector coefficient at @a ip.
846 const IntegrationPoint &ip) override;
847
848 /** @brief Evaluate the vector coefficients at all of the locations in the
849 integration rule and write the vectors into the columns of matrix @a
850 M. */
852 const IntegrationRule &ir) override;
853
854 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
855 /// the quadrature points.
856 ///
857 /// This function uses the efficient QuadratureFunction::ProjectGridFunction
858 /// to fill the QuadratureFunction.
859 void Project(QuadratureFunction &qf) override;
860
862};
863
864/// Vector coefficient defined as the Gradient of a scalar GridFunction
866{
867protected:
869
870public:
871
872 /** @brief Construct the coefficient with a scalar grid function @a gf. The
873 grid function is not owned by the coefficient. */
875
876 ///Set the scalar grid function.
877 void SetGridFunction(const GridFunction *gf);
878
879 ///Get the scalar grid function.
880 const GridFunction * GetGridFunction() const { return GridFunc; }
881
882 /// Evaluate the gradient vector coefficient at @a ip.
884 const IntegrationPoint &ip) override;
885
886 /** @brief Evaluate the gradient vector coefficient at all of the locations
887 in the integration rule and write the vectors into columns of matrix @a
888 M. */
890 const IntegrationRule &ir) override;
891
893};
894
895/// Vector coefficient defined as the Curl of a vector GridFunction
897{
898protected:
900
901public:
902 /** @brief Construct the coefficient with a vector grid function @a gf. The
903 grid function is not owned by the coefficient. */
905
906 /// Set the vector grid function.
907 void SetGridFunction(const GridFunction *gf);
908
909 /// Get the vector grid function.
910 const GridFunction * GetGridFunction() const { return GridFunc; }
911
913 /// Evaluate the vector curl coefficient at @a ip.
915 const IntegrationPoint &ip) override;
916
918};
919
920/// Scalar coefficient defined as the Divergence of a vector GridFunction
922{
923protected:
925
926public:
927 /** @brief Construct the coefficient with a vector grid function @a gf. The
928 grid function is not owned by the coefficient. */
930
931 /// Set the vector grid function.
932 void SetGridFunction(const GridFunction *gf) { GridFunc = gf; }
933
934 /// Get the vector grid function.
935 const GridFunction * GetGridFunction() const { return GridFunc; }
936
937 /// Evaluate the scalar divergence coefficient at @a ip.
939 const IntegrationPoint &ip) override;
940
942};
943
944/** @brief Vector coefficient defined by a scalar DeltaCoefficient and a
945 constant vector direction.
946
947 WARNING this cannot be used as a normal coefficient. The usual Eval method
948 is disabled. */
950{
951protected:
954
955public:
956 /// Construct with a vector of dimension @a vdim_.
958 : VectorCoefficient(vdim_), dir(vdim_), d() { }
959
960 /** @brief Construct with a Vector object representing the direction and a
961 unit delta function centered at (0.0,0.0,0.0) */
963 : VectorCoefficient(dir_.Size()), dir(dir_), d() { }
964
965 /** @brief Construct with a Vector object representing the direction and a
966 delta function scaled by @a s and centered at (x,0.0,0.0) */
968 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,s) { }
969
970 /** @brief Construct with a Vector object representing the direction and a
971 delta function scaled by @a s and centered at (x,y,0.0) */
973 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,s) { }
974
975 /** @brief Construct with a Vector object representing the direction and a
976 delta function scaled by @a s and centered at (x,y,z) */
978 real_t s)
979 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,z,s) { }
980
981 /// Set the time for internally stored coefficients
982 void SetTime(real_t t) override;
983
984 /// Replace the associated DeltaCoefficient with a new DeltaCoefficient.
985 /** The new DeltaCoefficient cannot have a specified weight Coefficient, i.e.
986 DeltaCoefficient::Weight() should return NULL. */
987 void SetDeltaCoefficient(const DeltaCoefficient& d_) { d = d_; }
988
989 /// Return the associated scalar DeltaCoefficient.
991
992 void SetScale(real_t s) { d.SetScale(s); }
993 void SetDirection(const Vector& d_);
994
995 void SetDeltaCenter(const Vector& center) { d.SetDeltaCenter(center); }
996 void GetDeltaCenter(Vector& center) { d.GetDeltaCenter(center); }
997
998 /** @brief Return the specified direction vector multiplied by the value
999 returned by DeltaCoefficient::EvalDelta() of the associated scalar
1000 DeltaCoefficient. */
1001 virtual void EvalDelta(Vector &V, ElementTransformation &T,
1002 const IntegrationPoint &ip);
1003
1005 /** @brief A VectorDeltaFunction cannot be evaluated. Calling this method
1006 will cause an MFEM error, terminating the application. */
1008 const IntegrationPoint &ip) override
1009 { mfem_error("VectorDeltaCoefficient::Eval"); }
1011};
1012
1013/** @brief Derived vector coefficient that has the value of the parent vector
1014 where it is active and is zero otherwise. */
1016{
1017private:
1019 Array<int> active_attr;
1020
1021public:
1022 /** @brief Construct with a parent vector coefficient and an array of zeros
1023 and ones representing the attributes for which this coefficient should be
1024 active. */
1027 { c = &vc; attr.Copy(active_attr); }
1028
1029 /// Set the time for internally stored coefficients
1030 void SetTime(real_t t) override;
1031
1032 /// Evaluate the vector coefficient at @a ip.
1033 void Eval(Vector &V, ElementTransformation &T,
1034 const IntegrationPoint &ip) override;
1035
1036 /** @brief Evaluate the vector coefficient at all of the locations in the
1037 integration rule and write the vectors into the columns of matrix @a
1038 M. */
1040 const IntegrationRule &ir) override;
1041};
1042
1044
1045/// Base class for Matrix Coefficients that optionally depend on time and space.
1047{
1048protected:
1051 bool symmetric; // deprecated
1052
1053public:
1054 /// Construct a dim x dim matrix coefficient.
1055 explicit MatrixCoefficient(int dim, bool symm=false)
1056 { height = width = dim; time = 0.; symmetric = symm; }
1057
1058 /// Construct a h x w matrix coefficient.
1059 MatrixCoefficient(int h, int w, bool symm=false) :
1060 height(h), width(w), time(0.), symmetric(symm) { }
1061
1062 /// Set the time for time dependent coefficients
1063 virtual void SetTime(real_t t) { time = t; }
1064
1065 /// Get the time for time dependent coefficients
1066 real_t GetTime() { return time; }
1067
1068 /// Get the height of the matrix.
1069 int GetHeight() const { return height; }
1070
1071 /// Get the width of the matrix.
1072 int GetWidth() const { return width; }
1073
1074 /// For backward compatibility get the width of the matrix.
1075 int GetVDim() const { return width; }
1076
1077 /** @deprecated Use SymmetricMatrixCoefficient instead */
1078 bool IsSymmetric() const { return symmetric; }
1079
1080 /** @brief Evaluate the matrix coefficient in the element described by @a T
1081 at the point @a ip, storing the result in @a K. */
1082 /** @note When this method is called, the caller must make sure that the
1083 IntegrationPoint associated with @a T is the same as @a ip. This can be
1084 achieved by calling T.SetIntPoint(&ip). */
1086 const IntegrationPoint &ip) = 0;
1087
1088 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
1089 /// the quadrature points. The matrix will be transposed or not according to
1090 /// the boolean argument @a transpose.
1091 ///
1092 /// The @a vdim of the QuadratureFunction should be equal to the height times
1093 /// the width of the matrix.
1094 virtual void Project(QuadratureFunction &qf, bool transpose=false);
1095
1096 /// (DEPRECATED) Evaluate a symmetric matrix coefficient.
1097 /** @brief Evaluate the upper triangular entries of the matrix coefficient
1098 in the symmetric case, similarly to Eval. Matrix entry (i,j) is stored
1099 in K[j - i + os_i] for 0 <= i <= j < width, os_0 = 0,
1100 os_{i+1} = os_i + width - i. That is, K = {M(0,0), ..., M(0,w-1),
1101 M(1,1), ..., M(1,w-1), ..., M(w-1,w-1) with w = width.
1102 @deprecated Use Eval() instead. */
1104 const IntegrationPoint &ip)
1105 { mfem_error("MatrixCoefficient::EvalSymmetric"); }
1106
1108};
1109
1110
1111/// A matrix coefficient that is constant in space and time.
1113{
1114private:
1115 DenseMatrix mat;
1116public:
1117 ///Construct using matrix @a m for the constant.
1119 : MatrixCoefficient(m.Height(), m.Width()), mat(m) { }
1121 /// Evaluate the matrix coefficient at @a ip.
1123 const IntegrationPoint &ip) override { M = mat; }
1124 /// Return a reference to the constant matrix.
1125 const DenseMatrix& GetMatrix() { return mat; }
1126};
1127
1128
1129/** @brief A piecewise matrix-valued coefficient with the pieces keyed off the
1130 element attribute numbers.
1131
1132 A value of zero will be returned for any missing attribute numbers.
1133
1134 This object will not assume ownership of any MatrixCoefficient
1135 objects passed to it. Consequently, the caller must ensure that
1136 the individual MatrixCoefficient objects are not deleted while
1137 this PWMatrixCoefficient is still in use.
1138
1139 \note The keys may either be domain attribute numbers or boundary
1140 attribute numbers. If the PWMatrixCoefficient is used with a
1141 domain integrator the keys are assumed to be domain attribute
1142 numbers. Similarly, if the PWMatrixCoefficient is used with a
1143 boundary integrator the keys are assumed to be boundary attribute
1144 numbers.
1145*/
1147{
1148private:
1149 /** Internal data structure to store pointers to the appropriate
1150 coefficients for different regions of the mesh. The keys used
1151 in the map are the mesh attribute numbers (either element
1152 attribute or boundary element attribute depending upon
1153 context). The values returned for any missing attributes will
1154 be zero. The coefficient pointers may be NULL in which case a
1155 value of zero is returned.
1156
1157 The MatrixCoefficient objects contained in this map are NOT
1158 owned by this PWMatrixCoefficient object. This means that they
1159 will not be deleted when this object is deleted also the caller
1160 must ensure that the various MatrixCoefficient objects are not
1161 deleted while this PWMatrixCoefficient is still needed.
1162 */
1163 std::map<int, MatrixCoefficient*> pieces;
1164
1165 /** Convenience function to check for compatible array lengths,
1166 loop over the arrays, and add their attribute/MatrixCoefficient
1167 pairs to the internal data structure.
1168 */
1169 void InitMap(const Array<int> & attr,
1170 const Array<MatrixCoefficient*> & coefs);
1171
1172public:
1173
1174 /// Constructs a piecewise matrix coefficient of dimension dim by dim
1175 explicit PWMatrixCoefficient(int dim, bool symm = false)
1176 : MatrixCoefficient(dim, symm) {}
1177
1178 /// Constructs a piecewise matrix coefficient of dimension h by w
1179 explicit PWMatrixCoefficient(int h, int w, bool symm = false)
1180 : MatrixCoefficient(h, w, symm) {}
1181
1182 /// Construct the coefficient using arrays describing the pieces
1183 /** \param dim - size of the square matrix-valued result
1184 \param attr - an array of attribute numbers for each piece
1185 \param coefs - the corresponding array of MatrixCoefficient pointers
1186 \param symm - true if the result will be symmetric, false otherwise
1187 Any missing attributes or NULL coefficient pointers will result in a
1188 zero matrix being returned.
1189
1190 \note Ownership of the MatrixCoefficient objects will NOT be
1191 transferred to this object.
1192 */
1194 const Array<MatrixCoefficient*> & coefs,
1195 bool symm=false)
1196 : MatrixCoefficient(dim, symm) { InitMap(attr, coefs); }
1197
1198 /// Construct the coefficient using arrays describing the pieces
1199 /** \param h - height of the matrix-valued result
1200 \param w - width of the matrix-valued result
1201 \param attr - an array of attribute numbers for each piece
1202 \param coefs - the corresponding array of MatrixCoefficient pointers
1203 \param symm - true if the result will be symmetric, false otherwise
1204 Any missing attributes or NULL coefficient pointers will result in a
1205 zero matrix being returned for that attribute.
1206
1207 \note Ownership of the MatrixCoefficient objects will NOT be
1208 transferred to this object.
1209 */
1210 PWMatrixCoefficient(int h, int w, const Array<int> & attr,
1211 const Array<MatrixCoefficient*> & coefs,
1212 bool symm=false)
1213 : MatrixCoefficient(h, w, symm) { InitMap(attr, coefs); }
1214
1215 /// Set the time for time dependent coefficients
1216 void SetTime(real_t t) override;
1217
1218 /// Replace a set of coefficients
1220 const Array<MatrixCoefficient*> & coefs)
1221 { InitMap(attr, coefs); }
1222
1223 /// Replace a single coefficient for a particular attribute
1224 void UpdateCoefficient(int attr, MatrixCoefficient & coef);
1225
1226 /// Remove a single MatrixCoefficient for a particular attribute
1227 void ZeroCoefficient(int attr)
1228 { pieces.erase(attr); }
1229
1230 /// Evaluate the coefficient.
1232 const IntegrationPoint &ip) override;
1233};
1234
1235/** @brief A matrix coefficient with an optional scalar coefficient multiplier
1236 \a q. The matrix function can either be represented by a std function or
1237 a constant matrix provided when constructing this object. */
1239{
1240private:
1241 std::function<void(const Vector &, DenseMatrix &)> Function;
1242 std::function<void(const Vector &, Vector &)> SymmFunction; // deprecated
1243 std::function<void(const Vector &, real_t, DenseMatrix &)> TDFunction;
1244
1245 Coefficient *Q;
1246 DenseMatrix mat;
1247
1248public:
1249 /// Define a time-independent square matrix coefficient from a std function
1250 /** \param dim - the size of the matrix
1251 \param F - time-independent function
1252 \param q - optional scalar Coefficient to scale the matrix coefficient */
1254 std::function<void(const Vector &, DenseMatrix &)> F,
1255 Coefficient *q = nullptr)
1256 : MatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1257 { }
1258
1259 /// Define a constant matrix coefficient times a scalar Coefficient
1260 /** \param m - constant matrix
1261 \param q - optional scalar Coefficient to scale the matrix coefficient */
1263 : MatrixCoefficient(m.Height(), m.Width()), Q(&q), mat(m)
1264 { }
1265
1266 /** @brief Define a time-independent symmetric square matrix coefficient from
1267 a std function */
1268 /** \param dim - the size of the matrix
1269 \param SymmF - function used in EvalSymmetric
1270 \param q - optional scalar Coefficient to scale the matrix coefficient
1271 @deprecated Use another constructor without setting SymmFunction. */
1273 std::function<void(const Vector &, Vector &)> SymmF,
1274 Coefficient *q = NULL)
1275 : MatrixCoefficient(dim, true), SymmFunction(std::move(SymmF)), Q(q), mat(0)
1276 { }
1277
1278 /// Define a time-dependent square matrix coefficient from a std function
1279 /** \param dim - the size of the matrix
1280 \param TDF - time-dependent function
1281 \param q - optional scalar Coefficient to scale the matrix coefficient */
1283 std::function<void(const Vector &, real_t, DenseMatrix &)> TDF,
1284 Coefficient *q = nullptr)
1285 : MatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1286 { }
1287
1288 /// Set the time for internally stored coefficients
1289 void SetTime(real_t t) override;
1290
1291 /// Evaluate the matrix coefficient at @a ip.
1293 const IntegrationPoint &ip) override;
1294
1295 /// (DEPRECATED) Evaluate the symmetric matrix coefficient at @a ip.
1296 /** @deprecated Use Eval() instead. */
1298 const IntegrationPoint &ip) override;
1299
1301};
1302
1303
1304/** @brief Matrix coefficient defined by a matrix of scalar coefficients.
1305 Coefficients that are not set will evaluate to zero in the vector. The
1306 coefficient is stored as a flat Array with indexing (i,j) -> i*width+j. */
1308{
1309private:
1311 Array<bool> ownCoeff;
1312
1313public:
1314 /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
1315 actual coefficients still need to be added with Set(). */
1316 explicit MatrixArrayCoefficient (int dim);
1317
1318 /// Set the time for internally stored coefficients
1319 void SetTime(real_t t) override;
1320
1321 /// Get the coefficient located at (i,j) in the matrix.
1322 Coefficient* GetCoeff (int i, int j) { return Coeff[i*width+j]; }
1323
1324 /** @brief Set the coefficient located at (i,j) in the matrix. By default by
1325 default this will take ownership of the Coefficient passed in, but this
1326 can be overridden with the @a own parameter. */
1327 void Set(int i, int j, Coefficient * c, bool own=true);
1328
1329 /// Set ownership of the coefficient at (i,j) in the matrix
1330 void SetOwnership(int i, int j, bool own) { ownCoeff[i*width+j] = own; }
1331
1332 /// Get ownership of the coefficient at (i,j) in the matrix
1333 bool GetOwnership(int i, int j) const { return ownCoeff[i*width+j]; }
1334
1336
1337 /// Evaluate coefficient located at (i,j) in the matrix using integration
1338 /// point @a ip.
1340 { return Coeff[i*width+j] ? Coeff[i*width+j] -> Eval(T, ip, GetTime()) : 0.0; }
1341
1342 /// Evaluate the matrix coefficient @a ip.
1344 const IntegrationPoint &ip) override;
1345
1346 virtual ~MatrixArrayCoefficient();
1347};
1348
1349/** @brief Matrix coefficient defined row-wise by an array of vector
1350 coefficients. Rows that are not set will evaluate to zero. The
1351 matrix coefficient is stored as an array indexing the rows of
1352 the matrix. */
1354{
1355private:
1357 Array<bool> ownCoeff;
1358
1359public:
1360 /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
1361 actual coefficients still need to be added with Set(). */
1362 explicit MatrixArrayVectorCoefficient (int dim);
1363
1364 /// Set the time for internally stored coefficients
1365 void SetTime(real_t t) override;
1366
1367 /// Get the vector coefficient located at the i-th row of the matrix
1368 VectorCoefficient* GetCoeff (int i) { return Coeff[i]; }
1369
1370 /** @brief Set the coefficient located at the i-th row of the matrix.
1371 By this will take ownership of the Coefficient passed in, but this
1372 can be overridden with the @a own parameter. */
1373 void Set(int i, VectorCoefficient * c, bool own=true);
1374
1375 /// Set ownership of the i'th coefficient
1376 void SetOwnership(int i, bool own) { ownCoeff[i] = own; }
1377
1378 /// Get ownership of the i'th coefficient
1379 bool GetOwnership(int i) const { return ownCoeff[i]; }
1380
1382
1383 /// Evaluate coefficient located at the i-th row of the matrix using integration
1384 /// point @a ip.
1385 void Eval(int i, Vector &V, ElementTransformation &T,
1386 const IntegrationPoint &ip);
1387
1388 /// Evaluate the matrix coefficient @a ip.
1390 const IntegrationPoint &ip) override;
1391
1393};
1394
1395
1396/** @brief Derived matrix coefficient that has the value of the parent matrix
1397 coefficient where it is active and is zero otherwise. */
1399{
1400private:
1402 Array<int> active_attr;
1403
1404public:
1405 /** @brief Construct with a parent matrix coefficient and an array of zeros
1406 and ones representing the attributes for which this coefficient should be
1407 active. */
1410 { c = &mc; attr.Copy(active_attr); }
1411
1412 /// Set the time for internally stored coefficients
1413 void SetTime(real_t t) override;
1414
1415 /// Evaluate the matrix coefficient at @a ip.
1417 const IntegrationPoint &ip) override;
1418};
1419
1420/// Coefficients based on sums, products, or other functions of coefficients.
1421///@{
1422/** @brief Scalar coefficient defined as the linear combination of two scalar
1423 coefficients or a scalar and a scalar coefficient */
1425{
1426private:
1427 real_t aConst;
1428 Coefficient * a;
1429 Coefficient * b;
1430
1431 real_t alpha;
1432 real_t beta;
1433
1434public:
1435 /// Constructor with one coefficient. Result is alpha_ * A + beta_ * B
1437 real_t alpha_ = 1.0, real_t beta_ = 1.0)
1438 : aConst(A), a(NULL), b(&B), alpha(alpha_), beta(beta_) { }
1439
1440 /// Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
1442 real_t alpha_ = 1.0, real_t beta_ = 1.0)
1443 : aConst(0.0), a(&A), b(&B), alpha(alpha_), beta(beta_) { }
1444
1445 /// Set the time for internally stored coefficients
1446 void SetTime(real_t t) override;
1447
1448 /// Reset the first term in the linear combination as a constant
1449 void SetAConst(real_t A) { a = NULL; aConst = A; }
1450 /// Return the first term in the linear combination
1451 real_t GetAConst() const { return aConst; }
1452
1453 /// Reset the first term in the linear combination
1454 void SetACoef(Coefficient &A) { a = &A; }
1455 /// Return the first term in the linear combination
1456 Coefficient * GetACoef() const { return a; }
1457
1458 /// Reset the second term in the linear combination
1459 void SetBCoef(Coefficient &B) { b = &B; }
1460 /// Return the second term in the linear combination
1461 Coefficient * GetBCoef() const { return b; }
1462
1463 /// Reset the factor in front of the first term in the linear combination
1464 void SetAlpha(real_t alpha_) { alpha = alpha_; }
1465 /// Return the factor in front of the first term in the linear combination
1466 real_t GetAlpha() const { return alpha; }
1467
1468 /// Reset the factor in front of the second term in the linear combination
1469 void SetBeta(real_t beta_) { beta = beta_; }
1470 /// Return the factor in front of the second term in the linear combination
1471 real_t GetBeta() const { return beta; }
1472
1473 /// Evaluate the coefficient at @a ip.
1475 const IntegrationPoint &ip) override
1476 {
1477 return alpha * ((a == NULL ) ? aConst : a->Eval(T, ip) )
1478 + beta * b->Eval(T, ip);
1479 }
1480};
1481
1482
1483/// Base class for symmetric matrix coefficients that optionally depend on time and space.
1485{
1486protected:
1487
1488 /// Internal matrix used when evaluating this coefficient as a DenseMatrix.
1490public:
1491 /// Construct a dim x dim matrix coefficient.
1494
1495 /// Get the size of the matrix.
1496 int GetSize() const { return height; }
1497
1498 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
1499 /// the quadrature points.
1500 ///
1501 /// @note As opposed to MatrixCoefficient::Project, this function stores only
1502 /// the @a symmetric part of the matrix at each quadrature point.
1503 ///
1504 /// The @a vdim of the coefficient should be equal to height*(height+1)/2.
1505 virtual void ProjectSymmetric(QuadratureFunction &qf);
1506
1507 /** @brief Evaluate the matrix coefficient in the element described by @a T
1508 at the point @a ip, storing the result as a symmetric matrix @a K. */
1509 /** @note When this method is called, the caller must make sure that the
1510 IntegrationPoint associated with @a T is the same as @a ip. This can be
1511 achieved by calling T.SetIntPoint(&ip). */
1513 const IntegrationPoint &ip) = 0;
1514
1515 /** @brief Evaluate the matrix coefficient in the element described by @a T
1516 at the point @a ip, storing the result as a dense matrix @a K. */
1517 /** This function allows the use of SymmetricMatrixCoefficient in situations
1518 where the symmetry is not taken advantage of.
1519
1520 @note When this method is called, the caller must make sure that the
1521 IntegrationPoint associated with @a T is the same as @a ip. This can be
1522 achieved by calling T.SetIntPoint(&ip). */
1524 const IntegrationPoint &ip) override;
1525
1526
1527 /// @deprecated Return a reference to the internal matrix used when evaluating this coefficient as a DenseMatrix.
1528 MFEM_DEPRECATED const DenseSymmetricMatrix& GetMatrix() { return mat_aux; }
1529
1531};
1532
1533
1534/// A matrix coefficient that is constant in space and time.
1536{
1537private:
1539
1540public:
1541 ///Construct using matrix @a m for the constant.
1545 /// Evaluate the matrix coefficient at @a ip.
1547 const IntegrationPoint &ip) override { M = mat; }
1548
1549 /// Return a reference to the constant matrix.
1550 const DenseSymmetricMatrix& GetMatrix() { return mat; }
1551
1552};
1553
1554
1555/** @brief A matrix coefficient with an optional scalar coefficient multiplier
1556 \a q. The matrix function can either be represented by a std function or
1557 a constant matrix provided when constructing this object. */
1559{
1560private:
1561 std::function<void(const Vector &, DenseSymmetricMatrix &)> Function;
1562 std::function<void(const Vector &, real_t, DenseSymmetricMatrix &)> TDFunction;
1563
1564 Coefficient *Q;
1566
1567public:
1568 /// Define a time-independent symmetric matrix coefficient from a std function
1569 /** \param dim - the size of the matrix
1570 \param F - time-independent function
1571 \param q - optional scalar Coefficient to scale the matrix coefficient */
1573 std::function<void(const Vector &, DenseSymmetricMatrix &)> F,
1574 Coefficient *q = nullptr)
1575 : SymmetricMatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1576 { }
1577
1578 /// Define a constant matrix coefficient times a scalar Coefficient
1579 /** \param m - constant matrix
1580 \param q - optional scalar Coefficient to scale the matrix coefficient */
1582 Coefficient &q)
1583 : SymmetricMatrixCoefficient(m.Height()), Q(&q), mat(m)
1584 { }
1585
1586 /// Define a time-dependent square matrix coefficient from a std function
1587 /** \param dim - the size of the matrix
1588 \param TDF - time-dependent function
1589 \param q - optional scalar Coefficient to scale the matrix coefficient */
1591 std::function<void(const Vector &, real_t, DenseSymmetricMatrix &)> TDF,
1592 Coefficient *q = nullptr)
1593 : SymmetricMatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1594 { }
1595
1596 /// Set the time for internally stored coefficients
1597 void SetTime(real_t t) override;
1598
1600 /// Evaluate the matrix coefficient at @a ip.
1602 const IntegrationPoint &ip) override;
1603
1605};
1606
1607
1608/** @brief Scalar coefficient defined as the product of two scalar coefficients
1609 or a scalar and a scalar coefficient. */
1611{
1612private:
1613 real_t aConst;
1614 Coefficient * a;
1615 Coefficient * b;
1616
1617public:
1618 /// Constructor with one coefficient. Result is A * B.
1620 : aConst(A), a(NULL), b(&B) { }
1621
1622 /// Constructor with two coefficients. Result is A * B.
1624 : aConst(0.0), a(&A), b(&B) { }
1625
1626 /// Set the time for internally stored coefficients
1627 void SetTime(real_t t) override;
1628
1629 /// Reset the first term in the product as a constant
1630 void SetAConst(real_t A) { a = NULL; aConst = A; }
1631 /// Return the first term in the product
1632 real_t GetAConst() const { return aConst; }
1633
1634 /// Reset the first term in the product
1635 void SetACoef(Coefficient &A) { a = &A; }
1636 /// Return the first term in the product
1637 Coefficient * GetACoef() const { return a; }
1638
1639 /// Reset the second term in the product
1640 void SetBCoef(Coefficient &B) { b = &B; }
1641 /// Return the second term in the product
1642 Coefficient * GetBCoef() const { return b; }
1643
1644 /// Evaluate the coefficient at @a ip.
1646 const IntegrationPoint &ip) override
1647 { return ((a == NULL ) ? aConst : a->Eval(T, ip) ) * b->Eval(T, ip); }
1648};
1649
1650/** @brief Scalar coefficient defined as the ratio of two scalars where one or
1651 both scalars are scalar coefficients. */
1653{
1654private:
1655 real_t aConst;
1656 real_t bConst;
1657 Coefficient * a;
1658 Coefficient * b;
1659
1660public:
1661 /** Initialize a coefficient which returns A / B where @a A is a
1662 constant and @a B is a scalar coefficient */
1664 : aConst(A), bConst(1.0), a(NULL), b(&B) { }
1665 /** Initialize a coefficient which returns A / B where @a A and @a B are both
1666 scalar coefficients */
1668 : aConst(0.0), bConst(1.0), a(&A), b(&B) { }
1669 /** Initialize a coefficient which returns A / B where @a A is a
1670 scalar coefficient and @a B is a constant */
1672 : aConst(0.0), bConst(B), a(&A), b(NULL) { }
1673
1674 /// Set the time for internally stored coefficients
1675 void SetTime(real_t t) override;
1676
1677 /// Reset the numerator in the ratio as a constant
1678 void SetAConst(real_t A) { a = NULL; aConst = A; }
1679 /// Return the numerator of the ratio
1680 real_t GetAConst() const { return aConst; }
1681
1682 /// Reset the denominator in the ratio as a constant
1683 void SetBConst(real_t B) { b = NULL; bConst = B; }
1684 /// Return the denominator of the ratio
1685 real_t GetBConst() const { return bConst; }
1686
1687 /// Reset the numerator in the ratio
1688 void SetACoef(Coefficient &A) { a = &A; }
1689 /// Return the numerator of the ratio
1690 Coefficient * GetACoef() const { return a; }
1691
1692 /// Reset the denominator in the ratio
1693 void SetBCoef(Coefficient &B) { b = &B; }
1694 /// Return the denominator of the ratio
1695 Coefficient * GetBCoef() const { return b; }
1696
1697 /// Evaluate the coefficient
1699 const IntegrationPoint &ip) override
1700 {
1701 real_t den = (b == NULL ) ? bConst : b->Eval(T, ip);
1702 MFEM_ASSERT(den != 0.0, "Division by zero in RatioCoefficient");
1703 return ((a == NULL ) ? aConst : a->Eval(T, ip) ) / den;
1704 }
1705};
1706
1707/// Scalar coefficient defined as a scalar raised to a power
1709{
1710private:
1711 Coefficient * a;
1712
1713 real_t p;
1714
1715public:
1716 /// Construct with a coefficient and a constant power @a p_. Result is A^p.
1718 : a(&A), p(p_) { }
1719
1720 /// Set the time for internally stored coefficients
1721 void SetTime(real_t t) override;
1722
1723 /// Reset the base coefficient
1724 void SetACoef(Coefficient &A) { a = &A; }
1725 /// Return the base coefficient
1726 Coefficient * GetACoef() const { return a; }
1727
1728 /// Reset the exponent
1729 void SetExponent(real_t p_) { p = p_; }
1730 /// Return the exponent
1731 real_t GetExponent() const { return p; }
1732
1733 /// Evaluate the coefficient at @a ip.
1735 const IntegrationPoint &ip) override
1736 { return pow(a->Eval(T, ip), p); }
1737};
1738
1739
1740/// Scalar coefficient defined as the inner product of two vector coefficients
1742{
1743private:
1746
1747 mutable Vector va;
1748 mutable Vector vb;
1749public:
1750 /// Construct with the two vector coefficients. Result is $ A \cdot B $.
1752
1753 /// Set the time for internally stored coefficients
1754 void SetTime(real_t t) override;
1755
1756 /// Reset the first vector in the inner product
1757 void SetACoef(VectorCoefficient &A) { a = &A; }
1758 /// Return the first vector coefficient in the inner product
1759 VectorCoefficient * GetACoef() const { return a; }
1760
1761 /// Reset the second vector in the inner product
1762 void SetBCoef(VectorCoefficient &B) { b = &B; }
1763 /// Return the second vector coefficient in the inner product
1764 VectorCoefficient * GetBCoef() const { return b; }
1765
1766 /// Evaluate the coefficient at @a ip.
1768 const IntegrationPoint &ip) override;
1769};
1770
1771/// Scalar coefficient defined as a cross product of two vectors in the xy-plane.
1773{
1774private:
1777
1778 mutable Vector va;
1779 mutable Vector vb;
1780
1781public:
1782 /// Constructor with two vector coefficients. Result is $ A_x B_y - A_y * B_x; $.
1784
1785 /// Set the time for internally stored coefficients
1786 void SetTime(real_t t) override;
1787
1788 /// Reset the first vector in the product
1789 void SetACoef(VectorCoefficient &A) { a = &A; }
1790 /// Return the first vector of the product
1791 VectorCoefficient * GetACoef() const { return a; }
1792
1793 /// Reset the second vector in the product
1794 void SetBCoef(VectorCoefficient &B) { b = &B; }
1795 /// Return the second vector of the product
1796 VectorCoefficient * GetBCoef() const { return b; }
1797
1798 /// Evaluate the coefficient at @a ip.
1800 const IntegrationPoint &ip) override;
1801};
1802
1803/// Scalar coefficient defined as the determinant of a matrix coefficient
1805{
1806private:
1808
1809 mutable DenseMatrix ma;
1810
1811public:
1812 /// Construct with the matrix.
1814
1815 /// Set the time for internally stored coefficients
1816 void SetTime(real_t t) override;
1817
1818 /// Reset the matrix coefficient
1819 void SetACoef(MatrixCoefficient &A) { a = &A; }
1820 /// Return the matrix coefficient
1821 MatrixCoefficient * GetACoef() const { return a; }
1822
1823 /// Evaluate the determinant coefficient at @a ip.
1825 const IntegrationPoint &ip) override;
1826};
1827
1828/// Scalar coefficient defined as the trace of a matrix coefficient
1830{
1831private:
1833
1834 mutable DenseMatrix ma;
1835
1836public:
1837 /// Construct with the matrix.
1839
1840 /// Set the time for internally stored coefficients
1841 void SetTime(real_t t) override;
1842
1843 /// Reset the matrix coefficient
1844 void SetACoef(MatrixCoefficient &A) { a = &A; }
1845 /// Return the matrix coefficient
1846 MatrixCoefficient * GetACoef() const { return a; }
1847
1848 /// Evaluate the trace coefficient at @a ip.
1850 const IntegrationPoint &ip) override;
1851};
1852
1853/// Vector coefficient defined as the linear combination of two vectors
1855{
1856private:
1857 VectorCoefficient * ACoef;
1858 VectorCoefficient * BCoef;
1859
1860 Vector A;
1861 Vector B;
1862
1863 Coefficient * alphaCoef;
1864 Coefficient * betaCoef;
1865
1866 real_t alpha;
1867 real_t beta;
1868
1869 mutable Vector va;
1870
1871public:
1872 /** Constructor with no coefficients.
1873 To be used with the various "Set" methods */
1875
1876 /** Constructor with two vector coefficients.
1877 Result is alpha_ * A + beta_ * B */
1879 real_t alpha_ = 1.0, real_t beta_ = 1.0);
1880
1881 /** Constructor with scalar coefficients.
1882 Result is alpha_ * A_ + beta_ * B_ */
1884 Coefficient &alpha_, Coefficient &beta_);
1885
1886 /// Set the time for internally stored coefficients
1887 void SetTime(real_t t) override;
1888
1889 /// Reset the first vector coefficient
1890 void SetACoef(VectorCoefficient &A_) { ACoef = &A_; }
1891 /// Return the first vector coefficient
1892 VectorCoefficient * GetACoef() const { return ACoef; }
1893
1894 /// Reset the second vector coefficient
1895 void SetBCoef(VectorCoefficient &B_) { BCoef = &B_; }
1896 /// Return the second vector coefficient
1897 VectorCoefficient * GetBCoef() const { return BCoef; }
1898
1899 /// Reset the factor in front of the first vector coefficient
1900 void SetAlphaCoef(Coefficient &A_) { alphaCoef = &A_; }
1901 /// Return the factor in front of the first vector coefficient
1902 Coefficient * GetAlphaCoef() const { return alphaCoef; }
1903
1904 /// Reset the factor in front of the second vector coefficient
1905 void SetBetaCoef(Coefficient &B_) { betaCoef = &B_; }
1906 /// Return the factor in front of the second vector coefficient
1907 Coefficient * GetBetaCoef() const { return betaCoef; }
1908
1909 /// Reset the first vector as a constant
1910 void SetA(const Vector &A_) { A = A_; ACoef = NULL; }
1911 /// Return the first vector constant
1912 const Vector & GetA() const { return A; }
1913
1914 /// Reset the second vector as a constant
1915 void SetB(const Vector &B_) { B = B_; BCoef = NULL; }
1916 /// Return the second vector constant
1917 const Vector & GetB() const { return B; }
1918
1919 /// Reset the factor in front of the first vector coefficient as a constant
1920 void SetAlpha(real_t alpha_) { alpha = alpha_; alphaCoef = NULL; }
1921 /// Return the factor in front of the first vector coefficient
1922 real_t GetAlpha() const { return alpha; }
1923
1924 /// Reset the factor in front of the second vector coefficient as a constant
1925 void SetBeta(real_t beta_) { beta = beta_; betaCoef = NULL; }
1926 /// Return the factor in front of the second vector coefficient
1927 real_t GetBeta() const { return beta; }
1928
1929 /// Evaluate the coefficient at @a ip.
1930 void Eval(Vector &V, ElementTransformation &T,
1931 const IntegrationPoint &ip) override;
1933};
1934
1935/// Vector coefficient defined as a product of scalar and vector coefficients.
1937{
1938private:
1939 real_t aConst;
1940 Coefficient * a;
1942
1943public:
1944 /// Constructor with constant and vector coefficient. Result is A * B.
1946
1947 /// Constructor with two coefficients. Result is A * B.
1949
1950 /// Set the time for internally stored coefficients
1951 void SetTime(real_t t) override;
1952
1953 /// Reset the scalar factor as a constant
1954 void SetAConst(real_t A) { a = NULL; aConst = A; }
1955 /// Return the scalar factor
1956 real_t GetAConst() const { return aConst; }
1957
1958 /// Reset the scalar factor
1959 void SetACoef(Coefficient &A) { a = &A; }
1960 /// Return the scalar factor
1961 Coefficient * GetACoef() const { return a; }
1962
1963 /// Reset the vector factor
1964 void SetBCoef(VectorCoefficient &B) { b = &B; }
1965 /// Return the vector factor
1966 VectorCoefficient * GetBCoef() const { return b; }
1967
1968 /// Evaluate the coefficient at @a ip.
1969 void Eval(Vector &V, ElementTransformation &T,
1970 const IntegrationPoint &ip) override;
1972};
1973
1974/// Vector coefficient defined as a normalized vector field (returns v/|v|)
1976{
1977private:
1979
1980 real_t tol;
1981
1982public:
1983 /** @brief Return a vector normalized to a length of one
1984
1985 This class evaluates the vector coefficient @a A and, if |A| > @a tol,
1986 returns the normalized vector A / |A|. If |A| <= @a tol, the zero
1987 vector is returned.
1988 */
1990
1991 /// Set the time for internally stored coefficients
1992 void SetTime(real_t t) override;
1993
1994 /// Reset the vector coefficient
1995 void SetACoef(VectorCoefficient &A) { a = &A; }
1996 /// Return the vector coefficient
1997 VectorCoefficient * GetACoef() const { return a; }
1998
1999 /// Evaluate the coefficient at @a ip.
2000 void Eval(Vector &V, ElementTransformation &T,
2001 const IntegrationPoint &ip) override;
2003};
2004
2005/// Vector coefficient defined as a cross product of two vectors
2007{
2008private:
2011
2012 mutable Vector va;
2013 mutable Vector vb;
2014
2015public:
2016 /// Construct with the two coefficients. Result is A x B.
2018
2019 /// Set the time for internally stored coefficients
2020 void SetTime(real_t t) override;
2021
2022 /// Reset the first term in the product
2023 void SetACoef(VectorCoefficient &A) { a = &A; }
2024 /// Return the first term in the product
2025 VectorCoefficient * GetACoef() const { return a; }
2026
2027 /// Reset the second term in the product
2028 void SetBCoef(VectorCoefficient &B) { b = &B; }
2029 /// Return the second term in the product
2030 VectorCoefficient * GetBCoef() const { return b; }
2031
2032 /// Evaluate the coefficient at @a ip.
2033 void Eval(Vector &V, ElementTransformation &T,
2034 const IntegrationPoint &ip) override;
2036};
2037
2038/** @brief Vector coefficient defined as a product of a matrix coefficient and
2039 a vector coefficient. */
2041{
2042private:
2045
2046 mutable DenseMatrix ma;
2047 mutable Vector vb;
2048
2049public:
2050 /// Constructor with two coefficients. Result is A*B.
2052
2053 /// Set the time for internally stored coefficients
2054 void SetTime(real_t t) override;
2055
2056 /// Reset the matrix coefficient
2057 void SetACoef(MatrixCoefficient &A) { a = &A; }
2058 /// Return the matrix coefficient
2059 MatrixCoefficient * GetACoef() const { return a; }
2060
2061 /// Reset the vector coefficient
2062 void SetBCoef(VectorCoefficient &B) { b = &B; }
2063 /// Return the vector coefficient
2064 VectorCoefficient * GetBCoef() const { return b; }
2065
2066 /// Evaluate the vector coefficient at @a ip.
2067 void Eval(Vector &V, ElementTransformation &T,
2068 const IntegrationPoint &ip) override;
2070};
2071
2072/// Convenient alias for the MatrixVectorProductCoefficient
2074
2075/// Constant matrix coefficient defined as the identity of dimension d
2077{
2078private:
2079 int dim;
2080
2081public:
2082 /// Construct with the dimension of the square identity matrix.
2085
2086 /// Evaluate the matrix coefficient at @a ip.
2088 const IntegrationPoint &ip) override;
2089};
2090
2091/// Matrix coefficient defined as the linear combination of two matrices
2093{
2094private:
2097
2098 real_t alpha;
2099 real_t beta;
2100
2101 mutable DenseMatrix ma;
2102
2103public:
2104 /// Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
2106 real_t alpha_ = 1.0, real_t beta_ = 1.0);
2107
2108 /// Set the time for internally stored coefficients
2109 void SetTime(real_t t) override;
2110
2111 /// Reset the first matrix coefficient
2112 void SetACoef(MatrixCoefficient &A) { a = &A; }
2113 /// Return the first matrix coefficient
2114 MatrixCoefficient * GetACoef() const { return a; }
2115
2116 /// Reset the second matrix coefficient
2117 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2118 /// Return the second matrix coefficient
2119 MatrixCoefficient * GetBCoef() const { return b; }
2120
2121 /// Reset the factor in front of the first matrix coefficient
2122 void SetAlpha(real_t alpha_) { alpha = alpha_; }
2123 /// Return the factor in front of the first matrix coefficient
2124 real_t GetAlpha() const { return alpha; }
2125
2126 /// Reset the factor in front of the second matrix coefficient
2127 void SetBeta(real_t beta_) { beta = beta_; }
2128 /// Return the factor in front of the second matrix coefficient
2129 real_t GetBeta() const { return beta; }
2130
2131 /// Evaluate the matrix coefficient at @a ip.
2133 const IntegrationPoint &ip) override;
2134};
2135
2136/// Matrix coefficient defined as the product of two matrices
2138{
2139private:
2142
2143 mutable DenseMatrix ma;
2144 mutable DenseMatrix mb;
2145
2146public:
2147 /// Construct with the two coefficients. Result is A * B.
2149
2150 /// Reset the first matrix coefficient
2151 void SetACoef(MatrixCoefficient &A) { a = &A; }
2152 /// Return the first matrix coefficient
2153 MatrixCoefficient * GetACoef() const { return a; }
2154
2155 /// Reset the second matrix coefficient
2156 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2157 /// Return the second matrix coefficient
2158 MatrixCoefficient * GetBCoef() const { return b; }
2159
2160 /// Evaluate the matrix coefficient at @a ip.
2162 const IntegrationPoint &ip) override;
2163};
2164
2165/** @brief Matrix coefficient defined as a product of a scalar coefficient and a
2166 matrix coefficient.*/
2168{
2169private:
2170 real_t aConst;
2171 Coefficient * a;
2173
2174public:
2175 /// Constructor with one coefficient. Result is A*B.
2177
2178 /// Constructor with two coefficients. Result is A*B.
2180
2181 /// Set the time for internally stored coefficients
2182 void SetTime(real_t t) override;
2183
2184 /// Reset the scalar factor as a constant
2185 void SetAConst(real_t A) { a = NULL; aConst = A; }
2186 /// Return the scalar factor
2187 real_t GetAConst() const { return aConst; }
2188
2189 /// Reset the scalar factor
2190 void SetACoef(Coefficient &A) { a = &A; }
2191 /// Return the scalar factor
2192 Coefficient * GetACoef() const { return a; }
2193
2194 /// Reset the matrix factor
2195 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2196 /// Return the matrix factor
2197 MatrixCoefficient * GetBCoef() const { return b; }
2198
2199 /// Evaluate the matrix coefficient at @a ip.
2201 const IntegrationPoint &ip) override;
2202};
2203
2204/// Matrix coefficient defined as the transpose of a matrix coefficient
2206{
2207private:
2209
2210public:
2211 /// Construct with the matrix coefficient. Result is $ A^T $.
2213
2214 /// Set the time for internally stored coefficients
2215 void SetTime(real_t t) override;
2216
2217 /// Reset the matrix coefficient
2218 void SetACoef(MatrixCoefficient &A) { a = &A; }
2219 /// Return the matrix coefficient
2220 MatrixCoefficient * GetACoef() const { return a; }
2221
2222 /// Evaluate the matrix coefficient at @a ip.
2224 const IntegrationPoint &ip) override;
2225};
2226
2227/// Matrix coefficient defined as the inverse of a matrix coefficient.
2229{
2230private:
2232
2233public:
2234 /// Construct with the matrix coefficient. Result is $ A^{-1} $.
2236
2237 /// Set the time for internally stored coefficients
2238 void SetTime(real_t t) override;
2239
2240 /// Reset the matrix coefficient
2241 void SetACoef(MatrixCoefficient &A) { a = &A; }
2242 /// Return the matrix coefficient
2243 MatrixCoefficient * GetACoef() const { return a; }
2244
2245 /// Evaluate the matrix coefficient at @a ip.
2247 const IntegrationPoint &ip) override;
2248};
2249
2250/// Matrix coefficient defined as the exponential of a matrix coefficient.
2252{
2253private:
2255
2256public:
2257 /// Construct the matrix coefficient. Result is $ \exp(A) $.
2259
2260 /// Set the time for internally stored coefficients
2261 void SetTime(real_t t) override;
2262
2263 /// Reset the matrix coefficient
2264 void SetACoef(MatrixCoefficient &A) { a = &A; }
2265 /// Return the matrix coefficient
2266 MatrixCoefficient * GetACoef() const { return a; }
2267
2268 /// Evaluate the matrix coefficient at @a ip.
2270 const IntegrationPoint &ip) override;
2271};
2272
2273/// Matrix coefficient defined as the outer product of two vector coefficients.
2275{
2276private:
2279
2280 mutable Vector va;
2281 mutable Vector vb;
2282
2283public:
2284 /// Construct with two vector coefficients. Result is $ A B^T $.
2286
2287 /// Set the time for internally stored coefficients
2288 void SetTime(real_t t) override;
2289
2290 /// Reset the first vector in the outer product
2291 void SetACoef(VectorCoefficient &A) { a = &A; }
2292 /// Return the first vector coefficient in the outer product
2293 VectorCoefficient * GetACoef() const { return a; }
2294
2295 /// Reset the second vector in the outer product
2296 void SetBCoef(VectorCoefficient &B) { b = &B; }
2297 /// Return the second vector coefficient in the outer product
2298 VectorCoefficient * GetBCoef() const { return b; }
2299
2300 /// Evaluate the matrix coefficient at @a ip.
2302 const IntegrationPoint &ip) override;
2303};
2304
2305/** @brief Matrix coefficient defined as -a k x k x, for a vector k and scalar a
2306
2307 This coefficient returns $a * (|k|^2 I - k \otimes k)$, where I is
2308 the identity matrix and $\otimes$ indicates the outer product. This
2309 can be evaluated for vectors of any dimension but in three
2310 dimensions it corresponds to computing the cross product with k twice.
2311*/
2313{
2314private:
2315 real_t aConst;
2316 Coefficient * a;
2318
2319 mutable Vector vk;
2320
2321public:
2324
2325 /// Set the time for internally stored coefficients
2326 void SetTime(real_t t) override;
2327
2328 /// Reset the scalar factor as a constant
2329 void SetAConst(real_t A) { a = NULL; aConst = A; }
2330 /// Return the scalar factor
2331 real_t GetAConst() const { return aConst; }
2332
2333 /// Reset the scalar factor
2334 void SetACoef(Coefficient &A) { a = &A; }
2335 /// Return the scalar factor
2336 Coefficient * GetACoef() const { return a; }
2337
2338 /// Reset the vector factor
2339 void SetKCoef(VectorCoefficient &K) { k = &K; }
2340 /// Return the vector factor
2341 VectorCoefficient * GetKCoef() const { return k; }
2342
2343 /// Evaluate the matrix coefficient at @a ip.
2345 const IntegrationPoint &ip) override;
2346};
2347///@}
2348
2349/** @brief Vector quadrature function coefficient which requires that the
2350 quadrature rules used for this vector coefficient be the same as those that
2351 live within the supplied QuadratureFunction. */
2353{
2354private:
2355 const QuadratureFunction &QuadF; //do not own
2356 int index;
2357
2358public:
2359 /// Constructor with a quadrature function as input
2361
2362 /** Set the starting index within the QuadFunc that'll be used to project
2363 outwards as well as the corresponding length. The projected length should
2364 have the bounds of 1 <= length <= (length QuadFunc - index). */
2365 void SetComponent(int index_, int length_);
2366
2367 const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2368
2370 void Eval(Vector &V, ElementTransformation &T,
2371 const IntegrationPoint &ip) override;
2372
2373 void Project(QuadratureFunction &qf) override;
2374
2376};
2377
2378/** @brief Quadrature function coefficient which requires that the quadrature
2379 rules used for this coefficient be the same as those that live within the
2380 supplied QuadratureFunction. */
2382{
2383private:
2384 const QuadratureFunction &QuadF;
2385
2386public:
2387 /// Constructor with a quadrature function as input
2389
2390 const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2391
2392 real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override;
2393
2394 void Project(QuadratureFunction &qf) override;
2395
2397};
2398
2399/// Flags that determine what storage optimizations to use in CoefficientVector
2400enum class CoefficientStorage : int
2401{
2402 FULL = 0, ///< Store the coefficient as a full QuadratureFunction.
2403 CONSTANTS = 1 << 0, ///< Store constants using only @a vdim entries.
2404 SYMMETRIC = 1 << 1, ///< Store the triangular part of symmetric matrices.
2405 COMPRESSED = CONSTANTS | SYMMETRIC ///< Enable all above compressions.
2406};
2407
2412
2414{
2415 return int(a) & int(b);
2416}
2417
2418
2419/// @brief Class to represent a coefficient evaluated at quadrature points.
2420///
2421/// In the general case, a CoefficientVector is the same as a QuadratureFunction
2422/// with a coefficient projected onto it.
2423///
2424/// This class allows for some "compression" of the coefficient data, according
2425/// to the storage flags given by CoefficientStorage. For example, constant
2426/// coefficients can be stored using only @a vdim values, and symmetric matrices
2427/// can be stored using e.g. the upper triangular part of the matrix.
2429{
2430protected:
2431 CoefficientStorage storage; ///< Storage optimizations (see CoefficientStorage).
2432 int vdim; ///< Number of values per quadrature point.
2433 QuadratureSpaceBase &qs; ///< Associated QuadratureSpaceBase.
2434 QuadratureFunction *qf; ///< Internal QuadratureFunction (owned, may be NULL).
2435public:
2436 /// Create an empty CoefficientVector.
2439
2440 /// @brief Create a CoefficientVector from the given Coefficient and
2441 /// QuadratureSpaceBase.
2442 ///
2443 /// If @a coeff is NULL, it will be interpreted as a constant with value one.
2444 /// @sa CoefficientStorage for a description of @a storage_.
2447
2448 /// @brief Create a CoefficientVector from the given Coefficient and
2449 /// QuadratureSpaceBase.
2450 ///
2451 /// @sa CoefficientStorage for a description of @a storage_.
2454
2455 /// @brief Create a CoefficientVector from the given VectorCoefficient and
2456 /// QuadratureSpaceBase.
2457 ///
2458 /// @sa CoefficientStorage for a description of @a storage_.
2461
2462 /// @brief Create a CoefficientVector from the given MatrixCoefficient and
2463 /// QuadratureSpaceBase.
2464 ///
2465 /// @sa CoefficientStorage for a description of @a storage_.
2468
2469 /// @brief Evaluate the given Coefficient at the quadrature points defined by
2470 /// @ref qs.
2471 void Project(Coefficient &coeff);
2472
2473 /// @brief Evaluate the given VectorCoefficient at the quadrature points
2474 /// defined by @ref qs.
2475 ///
2476 /// @sa CoefficientVector for a description of the @a compress argument.
2477 void Project(VectorCoefficient &coeff);
2478
2479 /// @brief Evaluate the given MatrixCoefficient at the quadrature points
2480 /// defined by @ref qs.
2481 ///
2482 /// @sa CoefficientVector for a description of the @a compress argument.
2483 void Project(MatrixCoefficient &coeff, bool transpose=false);
2484
2485 /// @brief Project the transpose of @a coeff.
2486 ///
2487 /// @sa Project(MatrixCoefficient&, QuadratureSpace&, bool, bool)
2489
2490 /// Make this vector a reference to the given QuadratureFunction.
2491 void MakeRef(const QuadratureFunction &qf_);
2492
2493 /// Set this vector to the given constant.
2494 void SetConstant(real_t constant);
2495
2496 /// Set this vector to the given constant vector.
2497 void SetConstant(const Vector &constant);
2498
2499 /// Set this vector to the given constant matrix.
2500 void SetConstant(const DenseMatrix &constant);
2501
2502 /// Set this vector to the given constant symmetric matrix.
2503 void SetConstant(const DenseSymmetricMatrix &constant);
2504
2505 /// Return the number of values per quadrature point.
2506 int GetVDim() const;
2507
2509};
2510
2511/** @brief Compute the Lp norm of a function f.
2512 $ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} $ */
2514 const IntegrationRule *irs[]);
2515
2516/** @brief Compute the Lp norm of a vector function f = {f_i}_i=1...N.
2517 $ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} $ */
2519 const IntegrationRule *irs[]);
2520
2521#ifdef MFEM_USE_MPI
2522/** @brief Compute the global Lp norm of a function f.
2523 $ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} $ */
2525 const IntegrationRule *irs[]);
2526
2527/** @brief Compute the global Lp norm of a vector function f = {f_i}_i=1...N.
2528 $ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} $ */
2530 const IntegrationRule *irs[]);
2531#endif
2532
2533}
2534
2535#endif
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:935
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
CartesianCoefficient(int comp_)
comp_ index of the desired component (0 -> x, 1 -> y, 2 -> z)
Scalar coefficient which returns the x-component of the evaluation point.
Scalar coefficient which returns the y-component of the evaluation point.
Scalar coefficient which returns the z-component of the evaluation point.
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 ~Coefficient()
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip, real_t t)
Evaluate the coefficient in the element described by T at the point ip at time t.
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.
ConstantCoefficient(real_t c=1.0)
c is value of constant function
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
Matrix coefficient defined as -a k x k x, for a vector k and scalar a.
void SetKCoef(VectorCoefficient &K)
Reset the vector factor.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
CrossCrossCoefficient(real_t A, VectorCoefficient &K)
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
real_t GetAConst() const
Return the scalar factor.
Coefficient * GetACoef() const
Return the scalar factor.
VectorCoefficient * GetKCoef() const
Return the vector factor.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
Vector coefficient defined as the Curl of a vector GridFunction.
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.
const GridFunction * GetGridFunction() const
Get the vector grid function.
CurlGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
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.
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
DeltaCoefficient(real_t x, real_t s)
Construct a delta function scaled by s and centered at (x,0.0,0.0)
const real_t * Center()
DeltaCoefficient(real_t x, real_t y, real_t s)
Construct a delta function scaled by s and centered at (x,y,0.0)
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
Coefficient * Weight()
See SetWeight() for description of the weight Coefficient.
real_t(* tdf)(real_t)
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTol(real_t tol_)
Set the tolerance used during projection onto GridFunction to identify the Mesh vertex where the Cent...
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error,...
real_t Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
DeltaCoefficient()
Construct a unit delta function centered at (0.0,0.0,0.0)
void SetDeltaCenter(const Vector &center)
Set the center location of the delta function.
DeltaCoefficient(real_t x, real_t y, real_t z, real_t s)
Construct a delta function scaled by s and centered at (x,y,z)
void SetFunction(real_t(*f)(real_t))
Set a time-dependent function that multiplies the Scale().
void SetWeight(Coefficient *w)
Set a weight Coefficient that multiplies the DeltaCoefficient.
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
void SetScale(real_t s_)
Set the scale value multiplying the delta function.
real_t Tol()
Return the tolerance used to identify the mesh vertices.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Scalar coefficient defined as the determinant of a matrix coefficient.
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the determinant coefficient at ip.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Scalar coefficient defined as the Divergence of a vector GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
const GridFunction * GetGridFunction() const
Get the vector grid function.
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...
Matrix coefficient defined as the exponential of a matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
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 .
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
A general function coefficient.
FunctionCoefficient(std::function< real_t(const Vector &, real_t)> TDF)
Define a time-dependent coefficient from a std function.
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
MFEM_DEPRECATED FunctionCoefficient(real_t(*f)(Vector &))
(DEPRECATED) Define a time-independent coefficient from a C-function
FunctionCoefficient(std::function< real_t(const Vector &)> F)
Define a time-independent coefficient from a std function.
MFEM_DEPRECATED FunctionCoefficient(real_t(*tdf)(Vector &, real_t))
(DEPRECATED) Define a time-dependent coefficient from a C-function
Vector coefficient defined as the Gradient of a scalar GridFunction.
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.
const GridFunction * GetGridFunction() const
Get 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...
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
const GridFunction * GetGridFunction() const
Get the internal GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the internal GridFunction.
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.
GridFunctionCoefficient(const GridFunction *gf, int comp=1)
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Constant matrix coefficient defined as the identity of dimension d.
IdentityMatrixCoefficient(int d)
Construct with the dimension of the square identity matrix.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
Scalar coefficient defined as the inner product of two vector coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the inner product.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the inner product.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the inner product.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the inner product.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Matrix coefficient defined as the inverse of a matrix coefficient.
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 .
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
Matrix coefficient defined by a matrix of scalar coefficients. Coefficients that are not set will eva...
Coefficient * GetCoeff(int i, int j)
Get the coefficient located at (i,j) in the matrix.
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...
bool GetOwnership(int i, int j) const
Get ownership of the coefficient at (i,j) in the matrix.
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 SetOwnership(int i, int j, bool own)
Set ownership of the coefficient at (i,j) in the matrix.
Matrix coefficient defined row-wise by an array of vector coefficients. Rows that are not set will ev...
void SetOwnership(int i, bool own)
Set ownership of the i'th coefficient.
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.
bool GetOwnership(int i) const
Get ownership of the i'th coefficient.
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...
VectorCoefficient * GetCoeff(int i)
Get the vector coefficient located at the i-th row of the matrix.
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....
MatrixCoefficient(int dim, bool symm=false)
Construct a dim x dim matrix coefficient.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t GetTime()
Get the time for time dependent coefficients.
int GetVDim() const
For backward compatibility get the width of the matrix.
MatrixCoefficient(int h, int w, bool symm=false)
Construct a h x w matrix coefficient.
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.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate a symmetric matrix coefficient.
int GetHeight() const
Get the height of the matrix.
A matrix coefficient that is constant in space and time.
MatrixConstantCoefficient(const DenseMatrix &m)
Construct using matrix m for the constant.
const DenseMatrix & GetMatrix()
Return a reference to the constant matrix.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent square matrix coefficient from a std function.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, DenseMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
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.
MatrixFunctionCoefficient(const DenseMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, Vector &)> SymmF, Coefficient *q=NULL)
Define a time-independent symmetric square matrix coefficient from a std function.
Matrix coefficient defined as the product of two matrices.
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
MatrixCoefficient * GetBCoef() const
Return the second matrix coefficient.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
Derived matrix coefficient that has the value of the parent matrix coefficient where it is active and...
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.
MatrixRestrictedCoefficient(MatrixCoefficient &mc, Array< int > &attr)
Construct with a parent matrix coefficient and an array of zeros and ones representing the attributes...
Matrix coefficient defined as the linear combination of two matrices.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixCoefficient * GetBCoef() const
Return the second matrix coefficient.
void SetBeta(real_t beta_)
Reset the factor in front of the second matrix coefficient.
real_t GetBeta() const
Return the factor in front of the second matrix coefficient.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first matrix coefficient.
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 SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
real_t GetAlpha() const
Return the factor in front of the first matrix coefficient.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
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.
void SetBCoef(VectorCoefficient &B)
Reset the vector coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
VectorCoefficient * GetBCoef() const
Return the vector coefficient.
Mesh data type.
Definition mesh.hpp:64
Vector coefficient defined as a normalized vector field (returns v/|v|)
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.
VectorCoefficient * GetACoef() const
Return the vector coefficient.
void SetACoef(VectorCoefficient &A)
Reset the vector coefficient.
NormalizedVectorCoefficient(VectorCoefficient &A, real_t tol=1e-6)
Return a vector normalized to a length of one.
Matrix coefficient defined as the outer product of two vector coefficients.
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the outer product.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the outer product.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the outer product.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the outer product.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
A piecewise coefficient with the pieces keyed off the element attribute numbers.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
PWCoefficient(const Array< int > &attr, const Array< Coefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
PWCoefficient()
Constructs a piecewise coefficient.
void ZeroCoefficient(int attr)
Remove a single Coefficient for a particular attribute.
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.
void UpdateCoefficients(const Array< int > &attr, const Array< Coefficient * > &coefs)
Replace a set of coefficients.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
int GetNConst()
Returns the number of constants representing different attributes.
real_t & operator()(int i)
Return a reference to the i-th constant.
void UpdateConstants(Vector &c)
Update the constants with vector c.
PWConstCoefficient(int NumOfSubD=0)
Constructs a piecewise constant coefficient in NumOfSubD subdomains.
PWConstCoefficient(Vector &c)
Construct the constant coefficient using a vector of constants.
void operator=(real_t c)
Set the constants for all attributes to constant c.
A piecewise matrix-valued coefficient with the pieces keyed off the element attribute numbers.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
PWMatrixCoefficient(int dim, bool symm=false)
Constructs a piecewise matrix coefficient of dimension dim by dim.
PWMatrixCoefficient(int dim, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
PWMatrixCoefficient(int h, int w, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
PWMatrixCoefficient(int h, int w, bool symm=false)
Constructs a piecewise matrix coefficient of dimension h by w.
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 ZeroCoefficient(int attr)
Remove a single MatrixCoefficient for a particular attribute.
void UpdateCoefficients(const Array< int > &attr, const Array< MatrixCoefficient * > &coefs)
Replace a set of coefficients.
A piecewise vector-valued coefficient with the pieces keyed off the element attribute numbers.
void UpdateCoefficients(const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Replace a set of coefficients.
PWVectorCoefficient(int vd, const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void ZeroCoefficient(int attr)
Remove a single VectorCoefficient for a particular attribute.
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.
PWVectorCoefficient(int vd)
Constructs a piecewise vector coefficient of dimension vd.
Class for parallel meshes.
Definition pmesh.hpp:34
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Scalar coefficient defined as a scalar raised to a power.
real_t GetExponent() const
Return the exponent.
void SetExponent(real_t p_)
Reset the exponent.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
Coefficient * GetACoef() const
Return the base coefficient.
PowerCoefficient(Coefficient &A, real_t p_)
Construct with a coefficient and a constant power p_. Result is A^p.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetACoef(Coefficient &A)
Reset the base coefficient.
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
void SetBCoef(Coefficient &B)
Reset the second term in the product.
Coefficient * GetACoef() const
Return the first term in the product.
ProductCoefficient(real_t A, Coefficient &B)
Constructor with one coefficient. Result is A * B.
ProductCoefficient(Coefficient &A, Coefficient &B)
Constructor with two coefficients. Result is A * B.
real_t GetAConst() const
Return the first term in the product.
Coefficient * GetBCoef() const
Return the second term in the product.
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 SetACoef(Coefficient &A)
Reset the first term in the product.
void SetAConst(real_t A)
Reset the first term in the product as a constant.
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.
const QuadratureFunction & GetQuadFunction() const
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition qspace.hpp:26
Scalar coefficient defined as the ratio of two scalars where one or both scalars are scalar coefficie...
Coefficient * GetBCoef() const
Return the denominator of the ratio.
void SetBConst(real_t B)
Reset the denominator in the ratio as a constant.
void SetAConst(real_t A)
Reset the numerator in the ratio as a constant.
void SetBCoef(Coefficient &B)
Reset the denominator in the ratio.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
RatioCoefficient(real_t A, Coefficient &B)
real_t GetBConst() const
Return the denominator of the ratio.
Coefficient * GetACoef() const
Return the numerator of the ratio.
RatioCoefficient(Coefficient &A, Coefficient &B)
RatioCoefficient(Coefficient &A, real_t B)
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t GetAConst() const
Return the numerator of the ratio.
void SetACoef(Coefficient &A)
Reset the numerator in the ratio.
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
RestrictedCoefficient(Coefficient &c_, Array< int > &attr)
Construct with a parent coefficient and an array with ones marking the attributes on which this coeff...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Matrix coefficient defined as a product of a scalar coefficient and a matrix coefficient.
real_t GetAConst() const
Return the scalar factor.
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetBCoef(MatrixCoefficient &B)
Reset the matrix factor.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixCoefficient * GetBCoef() const
Return the matrix factor.
ScalarMatrixProductCoefficient(real_t A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
Coefficient * GetACoef() const
Return the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
Vector coefficient defined as a product of scalar and vector coefficients.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t GetAConst() const
Return the scalar factor.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Coefficient * GetACoef() const
Return the scalar factor.
void SetBCoef(VectorCoefficient &B)
Reset the vector factor.
VectorCoefficient * GetBCoef() const
Return the vector factor.
ScalarVectorProductCoefficient(real_t A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
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.
Scalar coefficient defined as the linear combination of two scalar coefficients or a scalar and a sca...
void SetBCoef(Coefficient &B)
Reset the second term in the linear combination.
Coefficient * GetACoef() const
Return the first term in the linear combination.
SumCoefficient(real_t A, Coefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Constructor with one coefficient. Result is alpha_ * A + beta_ * B.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first term in the linear combination.
SumCoefficient(Coefficient &A, Coefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
void SetACoef(Coefficient &A)
Reset the first term in the linear combination.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetBeta(real_t beta_)
Reset the factor in front of the second term in the linear combination.
real_t GetBeta() const
Return the factor in front of the second term in the linear combination.
real_t GetAlpha() const
Return the factor in front of the first term in the linear combination.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Coefficient * GetBCoef() const
Return the second term in the linear combination.
void SetAConst(real_t A)
Reset the first term in the linear combination as a constant.
real_t GetAConst() const
Return the first term in the linear combination.
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 ...
SymmetricMatrixCoefficient(int dimension)
Construct a dim x dim matrix coefficient.
MFEM_DEPRECATED const DenseSymmetricMatrix & GetMatrix()
DenseSymmetricMatrix mat_aux
Internal matrix used when evaluating this coefficient as a DenseMatrix.
int GetSize() const
Get the size of the matrix.
A matrix coefficient that is constant in space and time.
const DenseSymmetricMatrix & GetMatrix()
Return a reference to the constant matrix.
SymmetricMatrixConstantCoefficient(const DenseSymmetricMatrix &m)
Construct using matrix m for the constant.
void Eval(DenseSymmetricMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
SymmetricMatrixFunctionCoefficient(const DenseSymmetricMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseSymmetricMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent symmetric matrix coefficient from a std function.
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, DenseSymmetricMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
Scalar coefficient defined as the trace of a matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the trace coefficient at ip.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
TraceCoefficient(MatrixCoefficient &A)
Construct with the matrix.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
A coefficient that depends on 1 or 2 parent coefficients and a transformation rule represented by a C...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
TransformedCoefficient(Coefficient *q, std::function< real_t(real_t)> F)
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
TransformedCoefficient(Coefficient *q1, Coefficient *q2, std::function< real_t(real_t, real_t)> F)
Matrix coefficient defined as the transpose of a matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
bool GetOwnership(int i) const
Get ownership of the i'th coefficient.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
Coefficient ** GetCoeffs()
Returns the entire array of coefficients.
void SetOwnership(int i, bool own)
Set ownership of the i'th 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)
Coefficient * GetCoeff(int i)
Returns i'th coefficient.
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.
VectorCoefficient(int vd)
Initialize the VectorCoefficient with vector dimension vd.
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.
const Vector & GetVec() const
Return a reference to the constant vector in this class.
VectorConstantCoefficient(const Vector &v)
Construct the coefficient with constant vector v.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Vector coefficient defined as a cross product of two vectors.
void SetBCoef(VectorCoefficient &B)
Reset the second term in the product.
VectorCoefficient * GetACoef() const
Return the first term in the product.
VectorCoefficient * GetBCoef() const
Return the second term in the product.
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 SetACoef(VectorCoefficient &A)
Reset the first term in the product.
Vector coefficient defined by a scalar DeltaCoefficient and a constant vector direction.
void GetDeltaCenter(Vector &center)
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t y, real_t z, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
VectorDeltaCoefficient(const Vector &dir_)
Construct with a Vector object representing the direction and a unit delta function centered at (0....
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetDeltaCoefficient(const DeltaCoefficient &d_)
Replace the associated DeltaCoefficient with a new DeltaCoefficient.
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...
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t y, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
DeltaCoefficient & GetDeltaCoefficient()
Return the associated scalar DeltaCoefficient.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
A VectorDeltaFunction cannot be evaluated. Calling this method will cause an MFEM error,...
VectorDeltaCoefficient(int vdim_)
Construct with a vector of dimension vdim_.
void SetDeltaCenter(const Vector &center)
A general vector function coefficient.
VectorFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, Vector &)> TDF, Coefficient *q=nullptr)
Define a time-dependent vector coefficient from a std function.
VectorFunctionCoefficient(int dim, std::function< void(const Vector &, Vector &)> F, Coefficient *q=nullptr)
Define a time-independent vector coefficient from a std function.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Vector coefficient defined by a vector GridFunction.
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.
const GridFunction * GetGridFunction() const
Returns a pointer to the grid function in this Coefficient.
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 ...
const QuadratureFunction & GetQuadFunction() const
void SetComponent(int index_, int length_)
Derived vector coefficient that has the value of the parent vector where it is active and is zero oth...
VectorRestrictedCoefficient(VectorCoefficient &vc, Array< int > &attr)
Construct with a parent vector coefficient and an array of zeros and ones representing the attributes...
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.
Scalar coefficient defined as a cross product of two vectors in the xy-plane.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the product.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the product.
VectorCoefficient * GetBCoef() const
Return the second vector of the product.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
VectorCoefficient * GetACoef() const
Return the first vector of the product.
Vector coefficient defined as the linear combination of two vectors.
void SetB(const Vector &B_)
Reset the second vector as a constant.
real_t GetAlpha() const
Return the factor in front of the first vector coefficient.
const Vector & GetB() const
Return the second vector constant.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first vector coefficient as a constant.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
const Vector & GetA() const
Return the first vector constant.
real_t GetBeta() const
Return the factor in front of the second vector coefficient.
void SetA(const Vector &A_)
Reset the first vector as a constant.
void SetBeta(real_t beta_)
Reset the factor in front of the second vector coefficient as a constant.
Coefficient * GetAlphaCoef() const
Return the factor in front of the first vector coefficient.
void SetBCoef(VectorCoefficient &B_)
Reset the second vector coefficient.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetAlphaCoef(Coefficient &A_)
Reset the factor in front of the first vector coefficient.
Coefficient * GetBetaCoef() const
Return the factor in front of the second vector coefficient.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient.
void SetBetaCoef(Coefficient &B_)
Reset the factor in front of the second vector coefficient.
VectorCoefficient * GetACoef() const
Return the first vector coefficient.
void SetACoef(VectorCoefficient &A_)
Reset the first vector coefficient.
Vector data type.
Definition vector.hpp:82
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
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition hooke.cpp:45
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
int operator&(CoefficientStorage a, CoefficientStorage b)
VectorCoefficient DiagonalMatrixCoefficient
void mfem_error(const char *msg)
Definition error.cpp:154
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
CartesianZCoefficient CylindricalZCoefficient
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.
@ COMPRESSED
Enable all above compressions.
@ FULL
Store the coefficient as a full QuadratureFunction.
CoefficientStorage operator|(CoefficientStorage a, CoefficientStorage b)
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
MatrixVectorProductCoefficient MatVecCoefficient
Convenient alias for the MatrixVectorProductCoefficient.
real_t p(const Vector &x, real_t t)