MFEM v4.9.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 /// Fill the QuadratureFunction @a qf with the piecewise constant values.
137 void Project(QuadratureFunction &qf) override;
138};
139
140/** @brief A piecewise coefficient with the pieces keyed off the element
141 attribute numbers.
142
143 A value of zero will be returned for any missing attribute numbers.
144
145 This object will not assume ownership of any Coefficient objects
146 passed to it. Consequently, the caller must ensure that the
147 individual Coefficient objects are not deleted while this
148 PWCoefficient is still in use.
149
150 \note The keys may either be domain attribute numbers or boundary
151 attribute numbers. If the PWCoefficient is used with a domain
152 integrator the keys are assumed to be domain attribute
153 numbers. Similarly, if the PWCoefficient is used with a boundary
154 integrator the keys are assumed to be boundary attribute numbers.
155*/
157{
158private:
159 /** Internal data structure to store pointers to the appropriate
160 coefficients for different regions of the mesh. The keys used
161 in the map are the mesh attribute numbers (either element
162 attribute or boundary element attribute depending upon
163 context). The values returned for any missing attributes will
164 be zero. The coefficient pointers may be NULL in which case a
165 value of zero is returned.
166
167 The Coefficient objects contained in this map are NOT owned by
168 this PWCoefficient object. This means that they will not be
169 deleted when this object is deleted also the caller must ensure
170 that the various Coefficient objects are not deleted while this
171 PWCoefficient is still needed.
172 */
173 std::map<int, Coefficient*> pieces;
174
175 /** Convenience function to check for compatible array lengths,
176 loop over the arrays, and add their attribute/Coefficient pairs
177 to the internal data structure.
178 */
179 void InitMap(const Array<int> & attr,
180 const Array<Coefficient*> & coefs);
181
182public:
183
184 /// Constructs a piecewise coefficient
185 explicit PWCoefficient() {}
186
187 /// Construct the coefficient using arrays describing the pieces
188 /** \param attr - an array of attribute numbers for each piece
189 \param coefs - the corresponding array of Coefficient pointers
190 Any missing attributes or NULL coefficient pointers will result in a
191 value of zero being returned for that attribute.
192
193 \note Ownership of the Coefficient objects will NOT be
194 transferred to this object.
195 */
197 const Array<Coefficient*> & coefs)
198 { InitMap(attr, coefs); }
199
200 /// Set the time for time dependent coefficients
201 void SetTime(real_t t) override;
202
203 /// Replace a set of coefficients
205 const Array<Coefficient*> & coefs)
206 { InitMap(attr, coefs); }
207
208 /// Replace a single Coefficient for a particular attribute
209 void UpdateCoefficient(int attr, Coefficient & coef)
210 { pieces[attr] = &coef; }
211
212 /// Remove a single Coefficient for a particular attribute
213 void ZeroCoefficient(int attr)
214 { pieces.erase(attr); }
215
216 /// Evaluate the coefficient.
218 const IntegrationPoint &ip) override;
219};
220
221/// A general function coefficient
223{
224protected:
225 std::function<real_t(const Vector &)> Function;
226 std::function<real_t(const Vector &, real_t)> TDFunction;
227
228public:
229 /// Define a time-independent coefficient from a std function
230 /** \param F time-independent std::function */
231 FunctionCoefficient(std::function<real_t(const Vector &)> F)
232 : Function(std::move(F))
233 { }
234
235 /// Define a time-dependent coefficient from a std function
236 /** \param TDF time-dependent function */
237 FunctionCoefficient(std::function<real_t(const Vector &, real_t)> TDF)
238 : TDFunction(std::move(TDF))
239 { }
240
241 /// (DEPRECATED) Define a time-independent coefficient from a C-function
242 /** @deprecated Use the method where the C-function, @a f, uses a const
243 Vector argument instead of Vector. */
244 MFEM_DEPRECATED FunctionCoefficient(real_t (*f)(Vector &))
245 {
246 // Cast first to (void*) to suppress a warning from newer version of
247 // Clang when using -Wextra.
248 Function = reinterpret_cast<real_t(*)(const Vector&)>((void*)f);
249 TDFunction = NULL;
250 }
251
252 /// (DEPRECATED) Define a time-dependent coefficient from a C-function
253 /** @deprecated Use the method where the C-function, @a tdf, uses a const
254 Vector argument instead of Vector. */
255 MFEM_DEPRECATED FunctionCoefficient(real_t (*tdf)(Vector &, real_t))
256 {
257 Function = NULL;
258 // Cast first to (void*) to suppress a warning from newer version of
259 // Clang when using -Wextra.
260 TDFunction =
261 reinterpret_cast<real_t(*)(const Vector&,real_t)>((void*)tdf);
262 }
263
264 /// Evaluate the coefficient at @a ip.
266 const IntegrationPoint &ip) override;
267};
268
269/// A common base class for returning individual components of the domain's
270/// Cartesian coordinates.
272{
273protected:
274 int comp;
276
277 /// @a comp_ index of the desired component (0 -> x, 1 -> y, 2 -> z)
278 CartesianCoefficient(int comp_) : comp(comp_), transip(3) {}
279
280public:
281 /// Evaluate the coefficient at @a ip.
283 const IntegrationPoint &ip) override;
284};
285
286/// Scalar coefficient which returns the x-component of the evaluation point
292
293/// Scalar coefficient which returns the y-component of the evaluation point
299
300/// Scalar coefficient which returns the z-component of the evaluation point
306
307/// Scalar coefficient which returns the radial distance from the axis of
308/// the evaluation point in the cylindrical coordinate system
310{
311private:
312 mutable Vector transip;
313
314public:
316
317 /// Evaluate the coefficient at @a ip.
319 const IntegrationPoint &ip) override;
320};
321
322/// Scalar coefficient which returns the angular position or azimuth (often
323/// denoted by theta) of the evaluation point in the cylindrical coordinate
324/// system
326{
327private:
328 mutable Vector transip;
329
330public:
332
333 /// Evaluate the coefficient at @a ip.
335 const IntegrationPoint &ip) override;
336};
337
338/// Scalar coefficient which returns the height or altitude of
339/// the evaluation point in the cylindrical coordinate system
341
342/// Scalar coefficient which returns the radial distance from the origin of
343/// the evaluation point in the spherical coordinate system
345{
346private:
347 mutable Vector transip;
348
349public:
351
352 /// Evaluate the coefficient at @a ip.
354 const IntegrationPoint &ip) override;
355};
356
357/// Scalar coefficient which returns the azimuthal angle (often denoted by phi)
358/// of the evaluation point in the spherical coordinate system
360{
361private:
362 mutable Vector transip;
363
364public:
366
367 /// Evaluate the coefficient at @a ip.
369 const IntegrationPoint &ip) override;
370};
371
372/// Scalar coefficient which returns the polar angle (often denoted by theta)
373/// of the evaluation point in the spherical coordinate system
375{
376private:
377 mutable Vector transip;
378
379public:
381
382 /// Evaluate the coefficient at @a ip.
384 const IntegrationPoint &ip) override;
385};
386
387class GridFunction;
388
389/// Coefficient defined by a GridFunction. This coefficient is mesh dependent.
391{
392private:
393 const GridFunction *GridF;
394 int Component;
395
396public:
397 GridFunctionCoefficient() : GridF(NULL), Component(1) { }
398 /** Construct GridFunctionCoefficient from a given GridFunction, and
399 optionally specify a component to use if it is a vector GridFunction. */
400 GridFunctionCoefficient (const GridFunction *gf, int comp = 1)
401 { GridF = gf; Component = comp; }
402
403 /// Set the internal GridFunction
404 void SetGridFunction(const GridFunction *gf) { GridF = gf; }
405
406 /// Get the internal GridFunction
407 const GridFunction * GetGridFunction() const { return GridF; }
408
409 /// Evaluate the coefficient at @a ip.
411 const IntegrationPoint &ip) override;
412
413 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
414 /// the quadrature points.
415 ///
416 /// This function uses the efficient QuadratureFunction::ProjectGridFunction
417 /// to fill the QuadratureFunction.
418 void Project(QuadratureFunction &qf) override;
419};
420
421
422/** @brief A coefficient that depends on 1 or 2 parent coefficients and a
423 transformation rule represented by a C-function.
424
425 $ C(x,t) = T(Q1(x,t)) $ or $ C(x,t) = T(Q1(x,t), Q2(x,t)) $
426
427 where T is the transformation rule, and Q1/Q2 are the parent coefficients.*/
429{
430private:
431 Coefficient * Q1;
432 Coefficient * Q2;
433 std::function<real_t(real_t)> Transform1;
434 std::function<real_t(real_t, real_t)> Transform2;
435
436public:
438 : Q1(q), Transform1(std::move(F)) { Q2 = 0; Transform2 = 0; }
440 std::function<real_t(real_t, real_t)> F)
441 : Q1(q1), Q2(q2), Transform2(std::move(F)) { Transform1 = 0; }
442
443 /// Set the time for internally stored coefficients
444 void SetTime(real_t t) override;
445
446 /// Evaluate the coefficient at @a ip.
447 real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override;
448};
449
450/** @brief Delta function coefficient optionally multiplied by a weight
451 coefficient and a scaled time dependent C-function.
452
453 $ F(x,t) = w(x,t) s T(t) d(x - xc) $
454
455 where w is the optional weight coefficient, @a s is a scale factor
456 T is an optional time-dependent function and d is a delta function.
457
458 WARNING this cannot be used as a normal coefficient. The usual Eval
459 method is disabled. */
461{
462protected:
465 int sdim;
467
468public:
469
470 /// Construct a unit delta function centered at (0.0,0.0,0.0)
472 {
473 center[0] = center[1] = center[2] = 0.; scale = 1.; tol = 1e-12;
474 weight = NULL; sdim = 0; tdf = NULL;
475 }
476
477 /// Construct a delta function scaled by @a s and centered at (x,0.0,0.0)
479 {
480 center[0] = x; center[1] = 0.; center[2] = 0.; scale = s; tol = 1e-12;
481 weight = NULL; sdim = 1; tdf = NULL;
482 }
483
484 /// Construct a delta function scaled by @a s and centered at (x,y,0.0)
486 {
487 center[0] = x; center[1] = y; center[2] = 0.; scale = s; tol = 1e-12;
488 weight = NULL; sdim = 2; tdf = NULL;
489 }
490
491 /// Construct a delta function scaled by @a s and centered at (x,y,z)
493 {
494 center[0] = x; center[1] = y; center[2] = z; scale = s; tol = 1e-12;
495 weight = NULL; sdim = 3; tdf = NULL;
496 }
497
498 /// Set the time for internally stored coefficients
499 void SetTime(real_t t) override;
500
501 /// Set the center location of the delta function.
502 void SetDeltaCenter(const Vector& center);
503
504 /// Set the scale value multiplying the delta function.
505 void SetScale(real_t s_) { scale = s_; }
506
507 /// Set a time-dependent function that multiplies the Scale().
508 void SetFunction(real_t (*f)(real_t)) { tdf = f; }
509
510 /** @brief Set the tolerance used during projection onto GridFunction to
511 identify the Mesh vertex where the Center() of the delta function
512 lies. (default 1e-12)*/
513 void SetTol(real_t tol_) { tol = tol_; }
514
515 /// Set a weight Coefficient that multiplies the DeltaCoefficient.
516 /** The weight Coefficient multiplies the value returned by EvalDelta() but
517 not the value returned by Scale().
518 The weight Coefficient is also used as the L2-weight function when
519 projecting the DeltaCoefficient onto a GridFunction, so that the weighted
520 integral of the projection is exactly equal to the Scale(). */
521 void SetWeight(Coefficient *w) { weight = w; }
522
523 /// Return a pointer to a c-array representing the center of the delta
524 /// function.
525 const real_t *Center() { return center; }
526
527 /** @brief Return the scale factor times the optional time dependent
528 function. Returns $ s T(t) $ with $ T(t) = 1 $ when
529 not set by the user. */
530 real_t Scale() { return tdf ? (*tdf)(GetTime())*scale : scale; }
531
532 /// Return the tolerance used to identify the mesh vertices
533 real_t Tol() { return tol; }
534
535 /// See SetWeight() for description of the weight Coefficient.
536 Coefficient *Weight() { return weight; }
537
538 /// Write the center of the delta function into @a center.
540
541 /// The value of the function assuming we are evaluating at the delta center.
543 /** @brief A DeltaFunction cannot be evaluated. Calling this method will
544 cause an MFEM error, terminating the application. */
546 { mfem_error("DeltaCoefficient::Eval"); return 0.; }
547 virtual ~DeltaCoefficient() { delete weight; }
548};
549
550/** @brief Derived coefficient that takes the value of the parent coefficient
551 for the active attributes and is zero otherwise. */
553{
554private:
555 Coefficient *c;
556 Array<int> active_attr;
557
558public:
559 /** @brief Construct with a parent coefficient and an array with
560 ones marking the attributes on which this coefficient should be
561 active. */
563 { c = &c_; attr.Copy(active_attr); }
564
565 /// Set the time for internally stored coefficients
566 void SetTime(real_t t) override;
567
568 /// Evaluate the coefficient at @a ip.
570 { return active_attr[T.Attribute-1] ? c->Eval(T, ip, GetTime()) : 0.0; }
571};
572
573/// Base class for vector Coefficients that optionally depend on time and space.
575{
576protected:
577 int vdim;
579
580public:
581 /// Initialize the VectorCoefficient with vector dimension @a vd.
582 VectorCoefficient(int vd) { vdim = vd; time = 0.; }
583
584 /// Set the time for time dependent coefficients
585 virtual void SetTime(real_t t) { time = t; }
586
587 /// Get the time for time dependent coefficients
588 real_t GetTime() { return time; }
589
590 /// Returns dimension of the vector.
591 int GetVDim() { return vdim; }
592
593 /** @brief Evaluate the vector coefficient in the element described by @a T
594 at the point @a ip, storing the result in @a V. */
595 /** @note When this method is called, the caller must make sure that the
596 IntegrationPoint associated with @a T is the same as @a ip. This can be
597 achieved by calling T.SetIntPoint(&ip). */
598 virtual void Eval(Vector &V, ElementTransformation &T,
599 const IntegrationPoint &ip) = 0;
600
601 /** @brief Evaluate the vector coefficient in the element described by @a T
602 at all points of @a ir, storing the result in @a M. */
603 /** The dimensions of @a M are GetVDim() by ir.GetNPoints() and they must be
604 set by the implementation of this method.
605
606 The general implementation provided by the base class (using the Eval
607 method for one IntegrationPoint at a time) can be overloaded for more
608 efficient implementation.
609
610 @note The IntegrationPoint associated with @a T is not used, and this
611 method will generally modify this IntegrationPoint associated with @a T.
612 */
613 virtual void Eval(DenseMatrix &M, ElementTransformation &T,
614 const IntegrationRule &ir);
615
616 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
617 /// the quadrature points.
618 ///
619 /// The @a vdim of the VectorCoefficient should be equal to the @a vdim of
620 /// the QuadratureFunction.
621 virtual void Project(QuadratureFunction &qf);
622
623 virtual ~VectorCoefficient() { }
624};
625
626
627/// Vector coefficient that is constant in space and time.
629{
630private:
631 Vector vec;
632public:
633 /// Construct the coefficient with constant vector @a v.
635 : VectorCoefficient(v.Size()), vec(v) { }
637
638 /// Evaluate the vector coefficient at @a ip.
640 const IntegrationPoint &ip) override { V = vec; }
641
642 /// Return a reference to the constant vector in this class.
643 const Vector& GetVec() const { return vec; }
644};
645
646/** @brief A piecewise vector-valued coefficient with the pieces keyed off the
647 element attribute numbers.
648
649 A value of zero will be returned for any missing attribute numbers.
650
651 This object will not assume ownership of any VectorCoefficient
652 objects passed to it. Consequently, the caller must ensure that
653 the individual VectorCoefficient objects are not deleted while
654 this PWVectorCoefficient is still in use.
655
656 \note The keys may either be domain attribute numbers or boundary
657 attribute numbers. If the PWVectorCoefficient is used with a
658 domain integrator the keys are assumed to be domain attribute
659 numbers. Similarly, if the PWVectorCoefficient is used with a
660 boundary integrator the keys are assumed to be boundary attribute
661 numbers.
662*/
664{
665private:
666 /** Internal data structure to store pointers to the appropriate
667 coefficients for different regions of the mesh. The keys used
668 in the map are the mesh attribute numbers (either element
669 attribute or boundary element attribute depending upon
670 context). The values returned for any missing attributes will
671 be zero. The coefficient pointers may be NULL in which case a
672 value of zero is returned.
673
674 The VectorCoefficient objects contained in this map are NOT
675 owned by this PWVectorCoefficient object. This means that they
676 will not be deleted when this object is deleted also the caller
677 must ensure that the various VectorCoefficient objects are not
678 deleted while this PWVectorCoefficient is still needed.
679 */
680 std::map<int, VectorCoefficient*> pieces;
681
682 /** Convenience function to check for compatible array lengths,
683 loop over the arrays, and add their attribute/VectorCoefficient
684 pairs to the internal data structure.
685 */
686 void InitMap(const Array<int> & attr,
687 const Array<VectorCoefficient*> & coefs);
688
689public:
690
691 /// Constructs a piecewise vector coefficient of dimension vd
692 explicit PWVectorCoefficient(int vd): VectorCoefficient(vd) {}
693
694 /// Construct the coefficient using arrays describing the pieces
695 /** \param vd - dimension of the vector-valued result
696 \param attr - an array of attribute numbers for each piece
697 \param coefs - the corresponding array of VectorCoefficient pointers
698 Any missing attributes or NULL coefficient pointers will result in a
699 zero vector being returned for that attribute.
700
701 \note Ownership of the VectorCoefficient objects will NOT be
702 transferred to this object.
703 */
704 PWVectorCoefficient(int vd, const Array<int> & attr,
705 const Array<VectorCoefficient*> & coefs)
706 : VectorCoefficient(vd) { InitMap(attr, coefs); }
707
708 /// Set the time for time dependent coefficients
709 void SetTime(real_t t) override;
710
711 /// Replace a set of coefficients
713 const Array<VectorCoefficient*> & coefs)
714 { InitMap(attr, coefs); }
715
716 /// Replace a single Coefficient for a particular attribute
717 void UpdateCoefficient(int attr, VectorCoefficient & coef);
718
719 /// Remove a single VectorCoefficient for a particular attribute
720 void ZeroCoefficient(int attr)
721 { pieces.erase(attr); }
722
723 /// Evaluate the coefficient.
725 const IntegrationPoint &ip) override;
727};
728
729/// A vector coefficient which returns the physical location of the
730/// evaluation point in the Cartesian coordinate system.
732{
733public:
734
736
738 /// Evaluate the vector coefficient at @a ip.
740 const IntegrationPoint &ip) override;
741
743};
744
745/// A general vector function coefficient
747{
748private:
749 std::function<void(const Vector &, Vector &)> Function;
750 std::function<void(const Vector &, real_t, Vector &)> TDFunction;
751 Coefficient *Q;
752
753public:
754 /// Define a time-independent vector coefficient from a std function
755 /** \param dim - the size of the vector
756 \param F - time-independent function
757 \param q - optional scalar Coefficient to scale the vector coefficient */
759 std::function<void(const Vector &, Vector &)> F,
760 Coefficient *q = nullptr)
761 : VectorCoefficient(dim), Function(std::move(F)), Q(q)
762 { }
763
764 /// Define a time-dependent vector coefficient from a std function
765 /** \param dim - the size of the vector
766 \param TDF - time-dependent function
767 \param q - optional scalar Coefficient to scale the vector coefficient */
769 std::function<void(const Vector &, real_t, Vector &)> TDF,
770 Coefficient *q = nullptr)
771 : VectorCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
772 { }
773
775 /// Evaluate the vector coefficient at @a ip.
777 const IntegrationPoint &ip) override;
778
780};
781
782/** @brief Vector coefficient defined by an array of scalar coefficients.
783 Coefficients that are not set will evaluate to zero in the vector. This
784 object takes ownership of the array of coefficients inside it and deletes
785 them at object destruction. */
787{
788private:
790 Array<bool> ownCoeff;
791
792public:
793 /** @brief Construct vector of dim coefficients. The actual coefficients
794 still need to be added with Set(). */
795 explicit VectorArrayCoefficient(int dim);
796
797 /// Set the time for internally stored coefficients
798 void SetTime(real_t t) override;
799
800 /// Returns i'th coefficient.
801 Coefficient* GetCoeff(int i) { return Coeff[i]; }
802
803 /// Returns the entire array of coefficients.
804 Coefficient **GetCoeffs() { return Coeff; }
805
806 /// Sets coefficient in the vector.
807 void Set(int i, Coefficient *c, bool own=true);
808
809 /// Set ownership of the i'th coefficient
810 void SetOwnership(int i, bool own) { ownCoeff[i] = own; }
811
812 /// Get ownership of the i'th coefficient
813 bool GetOwnership(int i) const { return ownCoeff[i]; }
814
815 /// Evaluates i'th component of the vector of coefficients and returns the
816 /// value.
818 { return Coeff[i] ? Coeff[i]->Eval(T, ip, GetTime()) : 0.0; }
819
821 /** @brief Evaluate the coefficient. Each element of vector V comes from the
822 associated array of scalar coefficients. */
824 const IntegrationPoint &ip) override;
825
826 /// Destroys vector coefficient.
827 virtual ~VectorArrayCoefficient();
828};
829
830/// Vector coefficient defined by a vector GridFunction
832{
833protected:
835
836public:
837 /** @brief Construct an empty coefficient. Calling Eval() before the grid
838 function is set will cause a segfault. */
840
841 /** @brief Construct the coefficient with grid function @a gf. The
842 grid function is not owned by the coefficient. */
844
845 /** @brief Set the grid function for this coefficient. Also sets the Vector
846 dimension to match that of the @a gf. */
847 void SetGridFunction(const GridFunction *gf);
848
849 /// Returns a pointer to the grid function in this Coefficient
850 const GridFunction * GetGridFunction() const { return GridFunc; }
851
852 /// Evaluate the vector coefficient at @a ip.
854 const IntegrationPoint &ip) override;
855
856 /** @brief Evaluate the vector coefficients at all of the locations in the
857 integration rule and write the vectors into the columns of matrix @a
858 M. */
860 const IntegrationRule &ir) override;
861
862 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
863 /// the quadrature points.
864 ///
865 /// This function uses the efficient QuadratureFunction::ProjectGridFunction
866 /// to fill the QuadratureFunction.
867 void Project(QuadratureFunction &qf) override;
868
870};
871
872/// Vector coefficient defined as the Gradient of a scalar GridFunction
874{
875protected:
877
878public:
879
880 /** @brief Construct the coefficient with a scalar grid function @a gf. The
881 grid function is not owned by the coefficient. */
883
884 ///Set the scalar grid function.
885 void SetGridFunction(const GridFunction *gf);
886
887 ///Get the scalar grid function.
888 const GridFunction * GetGridFunction() const { return GridFunc; }
889
890 /// Evaluate the gradient vector coefficient at @a ip.
892 const IntegrationPoint &ip) override;
893
894 /** @brief Evaluate the gradient vector coefficient at all of the locations
895 in the integration rule and write the vectors into columns of matrix @a
896 M. */
898 const IntegrationRule &ir) override;
899
900 /// @copydoc VectorCoefficient::Project(QuadratureFunction &)
901 void Project(QuadratureFunction &qf) override;
902
904};
905
906/// Vector coefficient defined as the Curl of a vector GridFunction
908{
909protected:
911
912public:
913 /** @brief Construct the coefficient with a vector grid function @a gf. The
914 grid function is not owned by the coefficient. */
916
917 /// Set the vector grid function.
918 void SetGridFunction(const GridFunction *gf);
919
920 /// Get the vector grid function.
921 const GridFunction * GetGridFunction() const { return GridFunc; }
922
924 /// Evaluate the vector curl coefficient at @a ip.
926 const IntegrationPoint &ip) override;
927
929};
930
931/// Scalar coefficient defined as the Divergence of a vector GridFunction
933{
934protected:
936
937public:
938 /** @brief Construct the coefficient with a vector grid function @a gf. The
939 grid function is not owned by the coefficient. */
941
942 /// Set the vector grid function.
943 void SetGridFunction(const GridFunction *gf) { GridFunc = gf; }
944
945 /// Get the vector grid function.
946 const GridFunction * GetGridFunction() const { return GridFunc; }
947
948 /// Evaluate the scalar divergence coefficient at @a ip.
950 const IntegrationPoint &ip) override;
951
953};
954
955/** @brief Vector coefficient defined by a scalar DeltaCoefficient and a
956 constant vector direction.
957
958 WARNING this cannot be used as a normal coefficient. The usual Eval method
959 is disabled. */
961{
962protected:
965
966public:
967 /// Construct with a vector of dimension @a vdim_.
969 : VectorCoefficient(vdim_), dir(vdim_), d() { }
970
971 /** @brief Construct with a Vector object representing the direction and a
972 unit delta function centered at (0.0,0.0,0.0) */
974 : VectorCoefficient(dir_.Size()), dir(dir_), d() { }
975
976 /** @brief Construct with a Vector object representing the direction and a
977 delta function scaled by @a s and centered at (x,0.0,0.0) */
979 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,s) { }
980
981 /** @brief Construct with a Vector object representing the direction and a
982 delta function scaled by @a s and centered at (x,y,0.0) */
984 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,s) { }
985
986 /** @brief Construct with a Vector object representing the direction and a
987 delta function scaled by @a s and centered at (x,y,z) */
989 real_t s)
990 : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,z,s) { }
991
992 /// Set the time for internally stored coefficients
993 void SetTime(real_t t) override;
994
995 /// Replace the associated DeltaCoefficient with a new DeltaCoefficient.
996 /** The new DeltaCoefficient cannot have a specified weight Coefficient, i.e.
997 DeltaCoefficient::Weight() should return NULL. */
998 void SetDeltaCoefficient(const DeltaCoefficient& d_) { d = d_; }
999
1000 /// Return the associated scalar DeltaCoefficient.
1002
1003 void SetScale(real_t s) { d.SetScale(s); }
1004 void SetDirection(const Vector& d_);
1005
1006 void SetDeltaCenter(const Vector& center) { d.SetDeltaCenter(center); }
1007 void GetDeltaCenter(Vector& center) { d.GetDeltaCenter(center); }
1008
1009 /** @brief Return the specified direction vector multiplied by the value
1010 returned by DeltaCoefficient::EvalDelta() of the associated scalar
1011 DeltaCoefficient. */
1012 virtual void EvalDelta(Vector &V, ElementTransformation &T,
1013 const IntegrationPoint &ip);
1014
1016 /** @brief A VectorDeltaFunction cannot be evaluated. Calling this method
1017 will cause an MFEM error, terminating the application. */
1019 const IntegrationPoint &ip) override
1020 { mfem_error("VectorDeltaCoefficient::Eval"); }
1022};
1023
1024/** @brief Derived vector coefficient that has the value of the parent vector
1025 where it is active and is zero otherwise. */
1027{
1028private:
1030 Array<int> active_attr;
1031
1032public:
1033 /** @brief Construct with a parent vector coefficient and an array of zeros
1034 and ones representing the attributes for which this coefficient should be
1035 active. */
1038 { c = &vc; attr.Copy(active_attr); }
1039
1040 /// Set the time for internally stored coefficients
1041 void SetTime(real_t t) override;
1042
1043 /// Evaluate the vector coefficient at @a ip.
1044 void Eval(Vector &V, ElementTransformation &T,
1045 const IntegrationPoint &ip) override;
1046
1047 /** @brief Evaluate the vector coefficient at all of the locations in the
1048 integration rule and write the vectors into the columns of matrix @a
1049 M. */
1051 const IntegrationRule &ir) override;
1052};
1053
1055
1056/// Base class for Matrix Coefficients that optionally depend on time and space.
1058{
1059protected:
1062 bool symmetric; // deprecated
1063
1064public:
1065 /// Construct a dim x dim matrix coefficient.
1066 explicit MatrixCoefficient(int dim, bool symm=false)
1067 { height = width = dim; time = 0.; symmetric = symm; }
1068
1069 /// Construct a h x w matrix coefficient.
1070 MatrixCoefficient(int h, int w, bool symm=false) :
1071 height(h), width(w), time(0.), symmetric(symm) { }
1072
1073 /// Set the time for time dependent coefficients
1074 virtual void SetTime(real_t t) { time = t; }
1075
1076 /// Get the time for time dependent coefficients
1077 real_t GetTime() { return time; }
1078
1079 /// Get the height of the matrix.
1080 int GetHeight() const { return height; }
1081
1082 /// Get the width of the matrix.
1083 int GetWidth() const { return width; }
1084
1085 /// For backward compatibility get the width of the matrix.
1086 int GetVDim() const { return width; }
1087
1088 /** @deprecated Use SymmetricMatrixCoefficient instead */
1089 bool IsSymmetric() const { return symmetric; }
1090
1091 /** @brief Evaluate the matrix coefficient in the element described by @a T
1092 at the point @a ip, storing the result in @a K. */
1093 /** @note When this method is called, the caller must make sure that the
1094 IntegrationPoint associated with @a T is the same as @a ip. This can be
1095 achieved by calling T.SetIntPoint(&ip). */
1097 const IntegrationPoint &ip) = 0;
1098
1099 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
1100 /// the quadrature points. The matrix will be transposed or not according to
1101 /// the boolean argument @a transpose.
1102 ///
1103 /// The @a vdim of the QuadratureFunction should be equal to the height times
1104 /// the width of the matrix.
1105 virtual void Project(QuadratureFunction &qf, bool transpose=false);
1106
1107 /// (DEPRECATED) Evaluate a symmetric matrix coefficient.
1108 /** @brief Evaluate the upper triangular entries of the matrix coefficient
1109 in the symmetric case, similarly to Eval. Matrix entry (i,j) is stored
1110 in K[j - i + os_i] for 0 <= i <= j < width, os_0 = 0,
1111 os_{i+1} = os_i + width - i. That is, K = {M(0,0), ..., M(0,w-1),
1112 M(1,1), ..., M(1,w-1), ..., M(w-1,w-1) with w = width.
1113 @deprecated Use Eval() instead. */
1115 const IntegrationPoint &ip)
1116 { mfem_error("MatrixCoefficient::EvalSymmetric"); }
1117
1119};
1120
1121
1122/// A matrix coefficient that is constant in space and time.
1124{
1125private:
1126 DenseMatrix mat;
1127public:
1128 ///Construct using matrix @a m for the constant.
1130 : MatrixCoefficient(m.Height(), m.Width()), mat(m) { }
1132 /// Evaluate the matrix coefficient at @a ip.
1134 const IntegrationPoint &ip) override { M = mat; }
1135 /// Return a reference to the constant matrix.
1136 const DenseMatrix& GetMatrix() { return mat; }
1137};
1138
1139
1140/** @brief A piecewise matrix-valued coefficient with the pieces keyed off the
1141 element attribute numbers.
1142
1143 A value of zero will be returned for any missing attribute numbers.
1144
1145 This object will not assume ownership of any MatrixCoefficient
1146 objects passed to it. Consequently, the caller must ensure that
1147 the individual MatrixCoefficient objects are not deleted while
1148 this PWMatrixCoefficient is still in use.
1149
1150 \note The keys may either be domain attribute numbers or boundary
1151 attribute numbers. If the PWMatrixCoefficient is used with a
1152 domain integrator the keys are assumed to be domain attribute
1153 numbers. Similarly, if the PWMatrixCoefficient is used with a
1154 boundary integrator the keys are assumed to be boundary attribute
1155 numbers.
1156*/
1158{
1159private:
1160 /** Internal data structure to store pointers to the appropriate
1161 coefficients for different regions of the mesh. The keys used
1162 in the map are the mesh attribute numbers (either element
1163 attribute or boundary element attribute depending upon
1164 context). The values returned for any missing attributes will
1165 be zero. The coefficient pointers may be NULL in which case a
1166 value of zero is returned.
1167
1168 The MatrixCoefficient objects contained in this map are NOT
1169 owned by this PWMatrixCoefficient object. This means that they
1170 will not be deleted when this object is deleted also the caller
1171 must ensure that the various MatrixCoefficient objects are not
1172 deleted while this PWMatrixCoefficient is still needed.
1173 */
1174 std::map<int, MatrixCoefficient*> pieces;
1175
1176 /** Convenience function to check for compatible array lengths,
1177 loop over the arrays, and add their attribute/MatrixCoefficient
1178 pairs to the internal data structure.
1179 */
1180 void InitMap(const Array<int> & attr,
1181 const Array<MatrixCoefficient*> & coefs);
1182
1183public:
1184
1185 /// Constructs a piecewise matrix coefficient of dimension dim by dim
1186 explicit PWMatrixCoefficient(int dim, bool symm = false)
1187 : MatrixCoefficient(dim, symm) {}
1188
1189 /// Constructs a piecewise matrix coefficient of dimension h by w
1190 explicit PWMatrixCoefficient(int h, int w, bool symm = false)
1191 : MatrixCoefficient(h, w, symm) {}
1192
1193 /// Construct the coefficient using arrays describing the pieces
1194 /** \param dim - size of the square matrix-valued result
1195 \param attr - an array of attribute numbers for each piece
1196 \param coefs - the corresponding array of MatrixCoefficient pointers
1197 \param symm - true if the result will be symmetric, false otherwise
1198 Any missing attributes or NULL coefficient pointers will result in a
1199 zero matrix being returned.
1200
1201 \note Ownership of the MatrixCoefficient objects will NOT be
1202 transferred to this object.
1203 */
1205 const Array<MatrixCoefficient*> & coefs,
1206 bool symm=false)
1207 : MatrixCoefficient(dim, symm) { InitMap(attr, coefs); }
1208
1209 /// Construct the coefficient using arrays describing the pieces
1210 /** \param h - height of the matrix-valued result
1211 \param w - width of the matrix-valued result
1212 \param attr - an array of attribute numbers for each piece
1213 \param coefs - the corresponding array of MatrixCoefficient pointers
1214 \param symm - true if the result will be symmetric, false otherwise
1215 Any missing attributes or NULL coefficient pointers will result in a
1216 zero matrix being returned for that attribute.
1217
1218 \note Ownership of the MatrixCoefficient objects will NOT be
1219 transferred to this object.
1220 */
1221 PWMatrixCoefficient(int h, int w, const Array<int> & attr,
1222 const Array<MatrixCoefficient*> & coefs,
1223 bool symm=false)
1224 : MatrixCoefficient(h, w, symm) { InitMap(attr, coefs); }
1225
1226 /// Set the time for time dependent coefficients
1227 void SetTime(real_t t) override;
1228
1229 /// Replace a set of coefficients
1231 const Array<MatrixCoefficient*> & coefs)
1232 { InitMap(attr, coefs); }
1233
1234 /// Replace a single coefficient for a particular attribute
1235 void UpdateCoefficient(int attr, MatrixCoefficient & coef);
1236
1237 /// Remove a single MatrixCoefficient for a particular attribute
1238 void ZeroCoefficient(int attr)
1239 { pieces.erase(attr); }
1240
1241 /// Evaluate the coefficient.
1243 const IntegrationPoint &ip) override;
1244};
1245
1246/** @brief A matrix coefficient with an optional scalar coefficient multiplier
1247 \a q. The matrix function can either be represented by a std function or
1248 a constant matrix provided when constructing this object. */
1250{
1251private:
1252 std::function<void(const Vector &, DenseMatrix &)> Function;
1253 std::function<void(const Vector &, Vector &)> SymmFunction; // deprecated
1254 std::function<void(const Vector &, real_t, DenseMatrix &)> TDFunction;
1255
1256 Coefficient *Q;
1257 DenseMatrix mat;
1258
1259public:
1260 /// Define a time-independent square matrix coefficient from a std function
1261 /** \param dim - the size of the matrix
1262 \param F - time-independent function
1263 \param q - optional scalar Coefficient to scale the matrix coefficient */
1265 std::function<void(const Vector &, DenseMatrix &)> F,
1266 Coefficient *q = nullptr)
1267 : MatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1268 { }
1269
1270 /// Define a constant matrix coefficient times a scalar Coefficient
1271 /** \param m - constant matrix
1272 \param q - optional scalar Coefficient to scale the matrix coefficient */
1274 : MatrixCoefficient(m.Height(), m.Width()), Q(&q), mat(m)
1275 { }
1276
1277 /** @brief Define a time-independent symmetric square matrix coefficient from
1278 a std function */
1279 /** \param dim - the size of the matrix
1280 \param SymmF - function used in EvalSymmetric
1281 \param q - optional scalar Coefficient to scale the matrix coefficient
1282 @deprecated Use another constructor without setting SymmFunction. */
1284 std::function<void(const Vector &, Vector &)> SymmF,
1285 Coefficient *q = NULL)
1286 : MatrixCoefficient(dim, true), SymmFunction(std::move(SymmF)), Q(q), mat(0)
1287 { }
1288
1289 /// Define a time-dependent square matrix coefficient from a std function
1290 /** \param dim - the size of the matrix
1291 \param TDF - time-dependent function
1292 \param q - optional scalar Coefficient to scale the matrix coefficient */
1294 std::function<void(const Vector &, real_t, DenseMatrix &)> TDF,
1295 Coefficient *q = nullptr)
1296 : MatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1297 { }
1298
1299 /// Set the time for internally stored coefficients
1300 void SetTime(real_t t) override;
1301
1302 /// Evaluate the matrix coefficient at @a ip.
1304 const IntegrationPoint &ip) override;
1305
1306 /// (DEPRECATED) Evaluate the symmetric matrix coefficient at @a ip.
1307 /** @deprecated Use Eval() instead. */
1309 const IntegrationPoint &ip) override;
1310
1312};
1313
1314
1315/** @brief Matrix coefficient defined by a matrix of scalar coefficients.
1316 Coefficients that are not set will evaluate to zero in the vector. The
1317 coefficient is stored as a flat Array with indexing (i,j) -> i*width+j. */
1319{
1320private:
1322 Array<bool> ownCoeff;
1323
1324public:
1325 /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
1326 actual coefficients still need to be added with Set(). */
1327 explicit MatrixArrayCoefficient (int dim);
1328
1329 /// Set the time for internally stored coefficients
1330 void SetTime(real_t t) override;
1331
1332 /// Get the coefficient located at (i,j) in the matrix.
1333 Coefficient* GetCoeff (int i, int j) { return Coeff[i*width+j]; }
1334
1335 /** @brief Set the coefficient located at (i,j) in the matrix. By default by
1336 default this will take ownership of the Coefficient passed in, but this
1337 can be overridden with the @a own parameter. */
1338 void Set(int i, int j, Coefficient * c, bool own=true);
1339
1340 /// Set ownership of the coefficient at (i,j) in the matrix
1341 void SetOwnership(int i, int j, bool own) { ownCoeff[i*width+j] = own; }
1342
1343 /// Get ownership of the coefficient at (i,j) in the matrix
1344 bool GetOwnership(int i, int j) const { return ownCoeff[i*width+j]; }
1345
1347
1348 /// Evaluate coefficient located at (i,j) in the matrix using integration
1349 /// point @a ip.
1351 { return Coeff[i*width+j] ? Coeff[i*width+j] -> Eval(T, ip, GetTime()) : 0.0; }
1352
1353 /// Evaluate the matrix coefficient @a ip.
1355 const IntegrationPoint &ip) override;
1356
1357 virtual ~MatrixArrayCoefficient();
1358};
1359
1360/** @brief Matrix coefficient defined row-wise by an array of vector
1361 coefficients. Rows that are not set will evaluate to zero. The
1362 matrix coefficient is stored as an array indexing the rows of
1363 the matrix. */
1365{
1366private:
1368 Array<bool> ownCoeff;
1369
1370public:
1371 /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
1372 actual coefficients still need to be added with Set(). */
1373 explicit MatrixArrayVectorCoefficient (int dim);
1374
1375 /// Set the time for internally stored coefficients
1376 void SetTime(real_t t) override;
1377
1378 /// Get the vector coefficient located at the i-th row of the matrix
1379 VectorCoefficient* GetCoeff (int i) { return Coeff[i]; }
1380
1381 /** @brief Set the coefficient located at the i-th row of the matrix.
1382 By this will take ownership of the Coefficient passed in, but this
1383 can be overridden with the @a own parameter. */
1384 void Set(int i, VectorCoefficient * c, bool own=true);
1385
1386 /// Set ownership of the i'th coefficient
1387 void SetOwnership(int i, bool own) { ownCoeff[i] = own; }
1388
1389 /// Get ownership of the i'th coefficient
1390 bool GetOwnership(int i) const { return ownCoeff[i]; }
1391
1393
1394 /// Evaluate coefficient located at the i-th row of the matrix using integration
1395 /// point @a ip.
1396 void Eval(int i, Vector &V, ElementTransformation &T,
1397 const IntegrationPoint &ip);
1398
1399 /// Evaluate the matrix coefficient @a ip.
1401 const IntegrationPoint &ip) override;
1402
1404};
1405
1406
1407/** @brief Derived matrix coefficient that has the value of the parent matrix
1408 coefficient where it is active and is zero otherwise. */
1410{
1411private:
1413 Array<int> active_attr;
1414
1415public:
1416 /** @brief Construct with a parent matrix coefficient and an array of zeros
1417 and ones representing the attributes for which this coefficient should be
1418 active. */
1421 { c = &mc; attr.Copy(active_attr); }
1422
1423 /// Set the time for internally stored coefficients
1424 void SetTime(real_t t) override;
1425
1426 /// Evaluate the matrix coefficient at @a ip.
1428 const IntegrationPoint &ip) override;
1429};
1430
1431/// Coefficients based on sums, products, or other functions of coefficients.
1432///@{
1433/** @brief Scalar coefficient defined as the linear combination of two scalar
1434 coefficients or a scalar and a scalar coefficient */
1436{
1437private:
1438 real_t aConst;
1439 Coefficient * a;
1440 Coefficient * b;
1441
1442 real_t alpha;
1443 real_t beta;
1444
1445public:
1446 /// Constructor with one coefficient. Result is alpha_ * A + beta_ * B
1448 real_t alpha_ = 1.0, real_t beta_ = 1.0)
1449 : aConst(A), a(NULL), b(&B), alpha(alpha_), beta(beta_) { }
1450
1451 /// Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
1453 real_t alpha_ = 1.0, real_t beta_ = 1.0)
1454 : aConst(0.0), a(&A), b(&B), alpha(alpha_), beta(beta_) { }
1455
1456 /// Set the time for internally stored coefficients
1457 void SetTime(real_t t) override;
1458
1459 /// @copydoc Coefficient::Project(QuadratureFunction &)
1460 void Project(QuadratureFunction &qf) override;
1461
1462 /// Reset the first term in the linear combination as a constant
1463 void SetAConst(real_t A) { a = NULL; aConst = A; }
1464 /// Return the first term in the linear combination
1465 real_t GetAConst() const { return aConst; }
1466
1467 /// Reset the first term in the linear combination
1468 void SetACoef(Coefficient &A) { a = &A; }
1469 /// Return the first term in the linear combination
1470 Coefficient * GetACoef() const { return a; }
1471
1472 /// Reset the second term in the linear combination
1473 void SetBCoef(Coefficient &B) { b = &B; }
1474 /// Return the second term in the linear combination
1475 Coefficient * GetBCoef() const { return b; }
1476
1477 /// Reset the factor in front of the first term in the linear combination
1478 void SetAlpha(real_t alpha_) { alpha = alpha_; }
1479 /// Return the factor in front of the first term in the linear combination
1480 real_t GetAlpha() const { return alpha; }
1481
1482 /// Reset the factor in front of the second term in the linear combination
1483 void SetBeta(real_t beta_) { beta = beta_; }
1484 /// Return the factor in front of the second term in the linear combination
1485 real_t GetBeta() const { return beta; }
1486
1487 /// Evaluate the coefficient at @a ip.
1489 const IntegrationPoint &ip) override
1490 {
1491 return alpha * ((a == NULL ) ? aConst : a->Eval(T, ip) )
1492 + beta * b->Eval(T, ip);
1493 }
1494};
1495
1496
1497/// Base class for symmetric matrix coefficients that optionally depend on time and space.
1499{
1500protected:
1501
1502 /// Internal matrix used when evaluating this coefficient as a DenseMatrix.
1504public:
1505 /// Construct a dim x dim matrix coefficient.
1508
1509 /// Get the size of the matrix.
1510 int GetSize() const { return height; }
1511
1512 /// @brief Fill the QuadratureFunction @a qf by evaluating the coefficient at
1513 /// the quadrature points.
1514 ///
1515 /// @note As opposed to MatrixCoefficient::Project, this function stores only
1516 /// the @a symmetric part of the matrix at each quadrature point.
1517 ///
1518 /// The @a vdim of the coefficient should be equal to height*(height+1)/2.
1519 virtual void ProjectSymmetric(QuadratureFunction &qf);
1520
1521 /** @brief Evaluate the matrix coefficient in the element described by @a T
1522 at the point @a ip, storing the result as a symmetric matrix @a K. */
1523 /** @note When this method is called, the caller must make sure that the
1524 IntegrationPoint associated with @a T is the same as @a ip. This can be
1525 achieved by calling T.SetIntPoint(&ip). */
1527 const IntegrationPoint &ip) = 0;
1528
1529 /** @brief Evaluate the matrix coefficient in the element described by @a T
1530 at the point @a ip, storing the result as a dense matrix @a K. */
1531 /** This function allows the use of SymmetricMatrixCoefficient in situations
1532 where the symmetry is not taken advantage of.
1533
1534 @note When this method is called, the caller must make sure that the
1535 IntegrationPoint associated with @a T is the same as @a ip. This can be
1536 achieved by calling T.SetIntPoint(&ip). */
1538 const IntegrationPoint &ip) override;
1539
1540
1541 /// @deprecated Return a reference to the internal matrix used when evaluating this coefficient as a DenseMatrix.
1542 MFEM_DEPRECATED const DenseSymmetricMatrix& GetMatrix() { return mat_aux; }
1543
1545};
1546
1547
1548/// A matrix coefficient that is constant in space and time.
1550{
1551private:
1553
1554public:
1555 ///Construct using matrix @a m for the constant.
1559 /// Evaluate the matrix coefficient at @a ip.
1561 const IntegrationPoint &ip) override { M = mat; }
1562
1563 /// Return a reference to the constant matrix.
1564 const DenseSymmetricMatrix& GetMatrix() { return mat; }
1565
1566};
1567
1568
1569/** @brief A matrix coefficient with an optional scalar coefficient multiplier
1570 \a q. The matrix function can either be represented by a std function or
1571 a constant matrix provided when constructing this object. */
1573{
1574private:
1575 std::function<void(const Vector &, DenseSymmetricMatrix &)> Function;
1576 std::function<void(const Vector &, real_t, DenseSymmetricMatrix &)> TDFunction;
1577
1578 Coefficient *Q;
1580
1581public:
1582 /// Define a time-independent symmetric matrix coefficient from a std function
1583 /** \param dim - the size of the matrix
1584 \param F - time-independent function
1585 \param q - optional scalar Coefficient to scale the matrix coefficient */
1587 std::function<void(const Vector &, DenseSymmetricMatrix &)> F,
1588 Coefficient *q = nullptr)
1589 : SymmetricMatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1590 { }
1591
1592 /// Define a constant matrix coefficient times a scalar Coefficient
1593 /** \param m - constant matrix
1594 \param q - optional scalar Coefficient to scale the matrix coefficient */
1596 Coefficient &q)
1597 : SymmetricMatrixCoefficient(m.Height()), Q(&q), mat(m)
1598 { }
1599
1600 /// Define a time-dependent square matrix coefficient from a std function
1601 /** \param dim - the size of the matrix
1602 \param TDF - time-dependent function
1603 \param q - optional scalar Coefficient to scale the matrix coefficient */
1605 std::function<void(const Vector &, real_t, DenseSymmetricMatrix &)> TDF,
1606 Coefficient *q = nullptr)
1607 : SymmetricMatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1608 { }
1609
1610 /// Set the time for internally stored coefficients
1611 void SetTime(real_t t) override;
1612
1614 /// Evaluate the matrix coefficient at @a ip.
1616 const IntegrationPoint &ip) override;
1617
1619};
1620
1621
1622/** @brief Scalar coefficient defined as the product of two scalar coefficients
1623 or a scalar and a scalar coefficient. */
1625{
1626private:
1627 real_t aConst;
1628 Coefficient * a;
1629 Coefficient * b;
1630
1631public:
1632 /// Constructor with one coefficient. Result is A * B.
1634 : aConst(A), a(NULL), b(&B) { }
1635
1636 /// Constructor with two coefficients. Result is A * B.
1638 : aConst(0.0), a(&A), b(&B) { }
1639
1640 /// Set the time for internally stored coefficients
1641 void SetTime(real_t t) override;
1642
1643 /// @copydoc Coefficient::Project(QuadratureFunction &)
1644 void Project(QuadratureFunction &qf) override;
1645
1646 /// Reset the first term in the product as a constant
1647 void SetAConst(real_t A) { a = NULL; aConst = A; }
1648 /// Return the first term in the product
1649 real_t GetAConst() const { return aConst; }
1650
1651 /// Reset the first term in the product
1652 void SetACoef(Coefficient &A) { a = &A; }
1653 /// Return the first term in the product
1654 Coefficient * GetACoef() const { return a; }
1655
1656 /// Reset the second term in the product
1657 void SetBCoef(Coefficient &B) { b = &B; }
1658 /// Return the second term in the product
1659 Coefficient * GetBCoef() const { return b; }
1660
1661 /// Evaluate the coefficient at @a ip.
1663 const IntegrationPoint &ip) override
1664 { return ((a == NULL ) ? aConst : a->Eval(T, ip) ) * b->Eval(T, ip); }
1665};
1666
1667/** @brief Scalar coefficient defined as the ratio of two scalars where one or
1668 both scalars are scalar coefficients. */
1670{
1671private:
1672 real_t aConst;
1673 real_t bConst;
1674 Coefficient * a;
1675 Coefficient * b;
1676
1677public:
1678 /** Initialize a coefficient which returns A / B where @a A is a
1679 constant and @a B is a scalar coefficient */
1681 : aConst(A), bConst(1.0), a(NULL), b(&B) { }
1682 /** Initialize a coefficient which returns A / B where @a A and @a B are both
1683 scalar coefficients */
1685 : aConst(0.0), bConst(1.0), a(&A), b(&B) { }
1686 /** Initialize a coefficient which returns A / B where @a A is a
1687 scalar coefficient and @a B is a constant */
1689 : aConst(0.0), bConst(B), a(&A), b(NULL) { }
1690
1691 /// Set the time for internally stored coefficients
1692 void SetTime(real_t t) override;
1693
1694 /// @copydoc Coefficient::Project(QuadratureFunction &)
1695 void Project(QuadratureFunction &qf) override;
1696
1697 /// Reset the numerator in the ratio as a constant
1698 void SetAConst(real_t A) { a = NULL; aConst = A; }
1699 /// Return the numerator of the ratio
1700 real_t GetAConst() const { return aConst; }
1701
1702 /// Reset the denominator in the ratio as a constant
1703 void SetBConst(real_t B) { b = NULL; bConst = B; }
1704 /// Return the denominator of the ratio
1705 real_t GetBConst() const { return bConst; }
1706
1707 /// Reset the numerator in the ratio
1708 void SetACoef(Coefficient &A) { a = &A; }
1709 /// Return the numerator of the ratio
1710 Coefficient * GetACoef() const { return a; }
1711
1712 /// Reset the denominator in the ratio
1713 void SetBCoef(Coefficient &B) { b = &B; }
1714 /// Return the denominator of the ratio
1715 Coefficient * GetBCoef() const { return b; }
1716
1717 /// Evaluate the coefficient
1719 const IntegrationPoint &ip) override
1720 {
1721 real_t den = (b == NULL ) ? bConst : b->Eval(T, ip);
1722 MFEM_ASSERT(den != 0.0, "Division by zero in RatioCoefficient");
1723 return ((a == NULL ) ? aConst : a->Eval(T, ip) ) / den;
1724 }
1725};
1726
1727/// Scalar coefficient defined as a scalar raised to a power
1729{
1730private:
1731 Coefficient * a;
1732
1733 real_t p;
1734
1735public:
1736 /// Construct with a coefficient and a constant power @a p_. Result is A^p.
1738 : a(&A), p(p_) { }
1739
1740 /// Set the time for internally stored coefficients
1741 void SetTime(real_t t) override;
1742
1743 /// Reset the base coefficient
1744 void SetACoef(Coefficient &A) { a = &A; }
1745 /// Return the base coefficient
1746 Coefficient * GetACoef() const { return a; }
1747
1748 /// Reset the exponent
1749 void SetExponent(real_t p_) { p = p_; }
1750 /// Return the exponent
1751 real_t GetExponent() const { return p; }
1752
1753 /// Evaluate the coefficient at @a ip.
1755 const IntegrationPoint &ip) override
1756 { return pow(a->Eval(T, ip), p); }
1757};
1758
1759
1760/// Scalar coefficient defined as the inner product of two vector coefficients
1762{
1763private:
1766
1767 mutable Vector va;
1768 mutable Vector vb;
1769public:
1770 /// Construct with the two vector coefficients. Result is $ A \cdot B $.
1772
1773 /// Set the time for internally stored coefficients
1774 void SetTime(real_t t) override;
1775
1776 /// Reset the first vector in the inner product
1777 void SetACoef(VectorCoefficient &A) { a = &A; }
1778 /// Return the first vector coefficient in the inner product
1779 VectorCoefficient * GetACoef() const { return a; }
1780
1781 /// Reset the second vector in the inner product
1782 void SetBCoef(VectorCoefficient &B) { b = &B; }
1783 /// Return the second vector coefficient in the inner product
1784 VectorCoefficient * GetBCoef() const { return b; }
1785
1786 /// Evaluate the coefficient at @a ip.
1788 const IntegrationPoint &ip) override;
1789
1790 /// @copydoc Coefficient::Project(QuadratureFunction &)
1791 void Project(QuadratureFunction &qf) override;
1792};
1793
1794/// Scalar coefficient defined as a cross product of two vectors in the xy-plane.
1796{
1797private:
1800
1801 mutable Vector va;
1802 mutable Vector vb;
1803
1804public:
1805 /// Constructor with two vector coefficients. Result is $ A_x B_y - A_y * B_x; $.
1807
1808 /// Set the time for internally stored coefficients
1809 void SetTime(real_t t) override;
1810
1811 /// Reset the first vector in the product
1812 void SetACoef(VectorCoefficient &A) { a = &A; }
1813 /// Return the first vector of the product
1814 VectorCoefficient * GetACoef() const { return a; }
1815
1816 /// Reset the second vector in the product
1817 void SetBCoef(VectorCoefficient &B) { b = &B; }
1818 /// Return the second vector of the product
1819 VectorCoefficient * GetBCoef() const { return b; }
1820
1821 /// Evaluate the coefficient at @a ip.
1823 const IntegrationPoint &ip) override;
1824};
1825
1826/// Scalar coefficient defined as the determinant of a matrix coefficient
1828{
1829private:
1831
1832 mutable DenseMatrix ma;
1833
1834public:
1835 /// Construct with the matrix.
1837
1838 /// Set the time for internally stored coefficients
1839 void SetTime(real_t t) override;
1840
1841 /// Reset the matrix coefficient
1842 void SetACoef(MatrixCoefficient &A) { a = &A; }
1843 /// Return the matrix coefficient
1844 MatrixCoefficient * GetACoef() const { return a; }
1845
1846 /// Evaluate the determinant coefficient at @a ip.
1848 const IntegrationPoint &ip) override;
1849};
1850
1851/// Scalar coefficient defined as the trace of a matrix coefficient
1853{
1854private:
1856
1857 mutable DenseMatrix ma;
1858
1859public:
1860 /// Construct with the matrix.
1862
1863 /// Set the time for internally stored coefficients
1864 void SetTime(real_t t) override;
1865
1866 /// Reset the matrix coefficient
1867 void SetACoef(MatrixCoefficient &A) { a = &A; }
1868 /// Return the matrix coefficient
1869 MatrixCoefficient * GetACoef() const { return a; }
1870
1871 /// Evaluate the trace coefficient at @a ip.
1873 const IntegrationPoint &ip) override;
1874};
1875
1876/// Vector coefficient defined as the linear combination of two vectors
1878{
1879private:
1880 VectorCoefficient * ACoef;
1881 VectorCoefficient * BCoef;
1882
1883 Vector A;
1884 Vector B;
1885
1886 Coefficient * alphaCoef;
1887 Coefficient * betaCoef;
1888
1889 real_t alpha;
1890 real_t beta;
1891
1892 mutable Vector va;
1893
1894public:
1895 /** Constructor with no coefficients.
1896 To be used with the various "Set" methods */
1898
1899 /** Constructor with two vector coefficients.
1900 Result is alpha_ * A + beta_ * B */
1902 real_t alpha_ = 1.0, real_t beta_ = 1.0);
1903
1904 /** Constructor with scalar coefficients.
1905 Result is alpha_ * A_ + beta_ * B_ */
1907 Coefficient &alpha_, Coefficient &beta_);
1908
1909 /// Set the time for internally stored coefficients
1910 void SetTime(real_t t) override;
1911
1912 /// Reset the first vector coefficient
1913 void SetACoef(VectorCoefficient &A_) { ACoef = &A_; }
1914 /// Return the first vector coefficient
1915 VectorCoefficient * GetACoef() const { return ACoef; }
1916
1917 /// Reset the second vector coefficient
1918 void SetBCoef(VectorCoefficient &B_) { BCoef = &B_; }
1919 /// Return the second vector coefficient
1920 VectorCoefficient * GetBCoef() const { return BCoef; }
1921
1922 /// Reset the factor in front of the first vector coefficient
1923 void SetAlphaCoef(Coefficient &A_) { alphaCoef = &A_; }
1924 /// Return the factor in front of the first vector coefficient
1925 Coefficient * GetAlphaCoef() const { return alphaCoef; }
1926
1927 /// Reset the factor in front of the second vector coefficient
1928 void SetBetaCoef(Coefficient &B_) { betaCoef = &B_; }
1929 /// Return the factor in front of the second vector coefficient
1930 Coefficient * GetBetaCoef() const { return betaCoef; }
1931
1932 /// Reset the first vector as a constant
1933 void SetA(const Vector &A_) { A = A_; ACoef = NULL; }
1934 /// Return the first vector constant
1935 const Vector & GetA() const { return A; }
1936
1937 /// Reset the second vector as a constant
1938 void SetB(const Vector &B_) { B = B_; BCoef = NULL; }
1939 /// Return the second vector constant
1940 const Vector & GetB() const { return B; }
1941
1942 /// Reset the factor in front of the first vector coefficient as a constant
1943 void SetAlpha(real_t alpha_) { alpha = alpha_; alphaCoef = NULL; }
1944 /// Return the factor in front of the first vector coefficient
1945 real_t GetAlpha() const { return alpha; }
1946
1947 /// Reset the factor in front of the second vector coefficient as a constant
1948 void SetBeta(real_t beta_) { beta = beta_; betaCoef = NULL; }
1949 /// Return the factor in front of the second vector coefficient
1950 real_t GetBeta() const { return beta; }
1951
1952 /// Evaluate the coefficient at @a ip.
1953 void Eval(Vector &V, ElementTransformation &T,
1954 const IntegrationPoint &ip) override;
1956};
1957
1958/// Vector coefficient defined as a product of scalar and vector coefficients.
1960{
1961private:
1962 real_t aConst;
1963 Coefficient * a;
1965
1966public:
1967 /// Constructor with constant and vector coefficient. Result is A * B.
1969
1970 /// Constructor with two coefficients. Result is A * B.
1972
1973 /// Set the time for internally stored coefficients
1974 void SetTime(real_t t) override;
1975
1976 /// Reset the scalar factor as a constant
1977 void SetAConst(real_t A) { a = NULL; aConst = A; }
1978 /// Return the scalar factor
1979 real_t GetAConst() const { return aConst; }
1980
1981 /// Reset the scalar factor
1982 void SetACoef(Coefficient &A) { a = &A; }
1983 /// Return the scalar factor
1984 Coefficient * GetACoef() const { return a; }
1985
1986 /// Reset the vector factor
1987 void SetBCoef(VectorCoefficient &B) { b = &B; }
1988 /// Return the vector factor
1989 VectorCoefficient * GetBCoef() const { return b; }
1990
1991 /// Evaluate the coefficient at @a ip.
1992 void Eval(Vector &V, ElementTransformation &T,
1993 const IntegrationPoint &ip) override;
1995};
1996
1997/// Vector coefficient defined as a normalized vector field (returns v/|v|)
1999{
2000private:
2002
2003 real_t tol;
2004
2005public:
2006 /** @brief Return a vector normalized to a length of one
2007
2008 This class evaluates the vector coefficient @a A and, if |A| > @a tol,
2009 returns the normalized vector A / |A|. If |A| <= @a tol, the zero
2010 vector is returned.
2011 */
2013
2014 /// Set the time for internally stored coefficients
2015 void SetTime(real_t t) override;
2016
2017 /// Reset the vector coefficient
2018 void SetACoef(VectorCoefficient &A) { a = &A; }
2019 /// Return the vector coefficient
2020 VectorCoefficient * GetACoef() const { return a; }
2021
2022 /// Evaluate the coefficient at @a ip.
2023 void Eval(Vector &V, ElementTransformation &T,
2024 const IntegrationPoint &ip) override;
2026};
2027
2028/// Vector coefficient defined as a cross product of two vectors
2030{
2031private:
2034
2035 mutable Vector va;
2036 mutable Vector vb;
2037
2038public:
2039 /// Construct with the two coefficients. Result is A x B.
2041
2042 /// Set the time for internally stored coefficients
2043 void SetTime(real_t t) override;
2044
2045 /// Reset the first term in the product
2046 void SetACoef(VectorCoefficient &A) { a = &A; }
2047 /// Return the first term in the product
2048 VectorCoefficient * GetACoef() const { return a; }
2049
2050 /// Reset the second term in the product
2051 void SetBCoef(VectorCoefficient &B) { b = &B; }
2052 /// Return the second term in the product
2053 VectorCoefficient * GetBCoef() const { return b; }
2054
2055 /// Evaluate the coefficient at @a ip.
2056 void Eval(Vector &V, ElementTransformation &T,
2057 const IntegrationPoint &ip) override;
2059};
2060
2061/** @brief Vector coefficient defined as a product of a matrix coefficient and
2062 a vector coefficient. */
2064{
2065private:
2068
2069 mutable DenseMatrix ma;
2070 mutable Vector vb;
2071
2072public:
2073 /// Constructor with two coefficients. Result is A*B.
2075
2076 /// Set the time for internally stored coefficients
2077 void SetTime(real_t t) override;
2078
2079 /// Reset the matrix coefficient
2080 void SetACoef(MatrixCoefficient &A) { a = &A; }
2081 /// Return the matrix coefficient
2082 MatrixCoefficient * GetACoef() const { return a; }
2083
2084 /// Reset the vector coefficient
2085 void SetBCoef(VectorCoefficient &B) { b = &B; }
2086 /// Return the vector coefficient
2087 VectorCoefficient * GetBCoef() const { return b; }
2088
2089 /// Evaluate the vector coefficient at @a ip.
2090 void Eval(Vector &V, ElementTransformation &T,
2091 const IntegrationPoint &ip) override;
2093};
2094
2095/// Convenient alias for the MatrixVectorProductCoefficient
2097
2098/// Constant matrix coefficient defined as the identity of dimension d
2100{
2101private:
2102 int dim;
2103
2104public:
2105 /// Construct with the dimension of the square identity matrix.
2108
2109 /// Evaluate the matrix coefficient at @a ip.
2111 const IntegrationPoint &ip) override;
2112};
2113
2114/// Matrix coefficient defined as the linear combination of two matrices
2116{
2117private:
2120
2121 real_t alpha;
2122 real_t beta;
2123
2124 mutable DenseMatrix ma;
2125
2126public:
2127 /// Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
2129 real_t alpha_ = 1.0, real_t beta_ = 1.0);
2130
2131 /// Set the time for internally stored coefficients
2132 void SetTime(real_t t) override;
2133
2134 /// Reset the first matrix coefficient
2135 void SetACoef(MatrixCoefficient &A) { a = &A; }
2136 /// Return the first matrix coefficient
2137 MatrixCoefficient * GetACoef() const { return a; }
2138
2139 /// Reset the second matrix coefficient
2140 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2141 /// Return the second matrix coefficient
2142 MatrixCoefficient * GetBCoef() const { return b; }
2143
2144 /// Reset the factor in front of the first matrix coefficient
2145 void SetAlpha(real_t alpha_) { alpha = alpha_; }
2146 /// Return the factor in front of the first matrix coefficient
2147 real_t GetAlpha() const { return alpha; }
2148
2149 /// Reset the factor in front of the second matrix coefficient
2150 void SetBeta(real_t beta_) { beta = beta_; }
2151 /// Return the factor in front of the second matrix coefficient
2152 real_t GetBeta() const { return beta; }
2153
2154 /// Evaluate the matrix coefficient at @a ip.
2156 const IntegrationPoint &ip) override;
2157};
2158
2159/// Matrix coefficient defined as the product of two matrices
2161{
2162private:
2165
2166 mutable DenseMatrix ma;
2167 mutable DenseMatrix mb;
2168
2169public:
2170 /// Construct with the two coefficients. Result is A * B.
2172
2173 /// Reset the first matrix coefficient
2174 void SetACoef(MatrixCoefficient &A) { a = &A; }
2175 /// Return the first matrix coefficient
2176 MatrixCoefficient * GetACoef() const { return a; }
2177
2178 /// Reset the second matrix coefficient
2179 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2180 /// Return the second matrix coefficient
2181 MatrixCoefficient * GetBCoef() const { return b; }
2182
2183 /// Evaluate the matrix coefficient at @a ip.
2185 const IntegrationPoint &ip) override;
2186};
2187
2188/** @brief Matrix coefficient defined as a product of a scalar coefficient and a
2189 matrix coefficient.*/
2191{
2192private:
2193 real_t aConst;
2194 Coefficient * a;
2196
2197public:
2198 /// Constructor with one coefficient. Result is A*B.
2200
2201 /// Constructor with two coefficients. Result is A*B.
2203
2204 /// Set the time for internally stored coefficients
2205 void SetTime(real_t t) override;
2206
2207 /// Reset the scalar factor as a constant
2208 void SetAConst(real_t A) { a = NULL; aConst = A; }
2209 /// Return the scalar factor
2210 real_t GetAConst() const { return aConst; }
2211
2212 /// Reset the scalar factor
2213 void SetACoef(Coefficient &A) { a = &A; }
2214 /// Return the scalar factor
2215 Coefficient * GetACoef() const { return a; }
2216
2217 /// Reset the matrix factor
2218 void SetBCoef(MatrixCoefficient &B) { b = &B; }
2219 /// Return the matrix factor
2220 MatrixCoefficient * GetBCoef() const { return b; }
2221
2222 /// Evaluate the matrix coefficient at @a ip.
2224 const IntegrationPoint &ip) override;
2225};
2226
2227/// Matrix coefficient defined as the transpose of a matrix coefficient
2229{
2230private:
2232
2233public:
2234 /// Construct with the matrix coefficient. Result is $ A^T $.
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 inverse of a matrix coefficient.
2252{
2253private:
2255
2256public:
2257 /// Construct with the matrix coefficient. Result is $ A^{-1} $.
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 exponential of a matrix coefficient.
2275{
2276private:
2278
2279public:
2280 /// Construct the matrix coefficient. Result is $ \exp(A) $.
2282
2283 /// Set the time for internally stored coefficients
2284 void SetTime(real_t t) override;
2285
2286 /// Reset the matrix coefficient
2287 void SetACoef(MatrixCoefficient &A) { a = &A; }
2288 /// Return the matrix coefficient
2289 MatrixCoefficient * GetACoef() const { return a; }
2290
2291 /// Evaluate the matrix coefficient at @a ip.
2293 const IntegrationPoint &ip) override;
2294};
2295
2296/// Matrix coefficient defined as the outer product of two vector coefficients.
2298{
2299private:
2302
2303 mutable Vector va;
2304 mutable Vector vb;
2305
2306public:
2307 /// Construct with two vector coefficients. Result is $ A B^T $.
2309
2310 /// Set the time for internally stored coefficients
2311 void SetTime(real_t t) override;
2312
2313 /// Reset the first vector in the outer product
2314 void SetACoef(VectorCoefficient &A) { a = &A; }
2315 /// Return the first vector coefficient in the outer product
2316 VectorCoefficient * GetACoef() const { return a; }
2317
2318 /// Reset the second vector in the outer product
2319 void SetBCoef(VectorCoefficient &B) { b = &B; }
2320 /// Return the second vector coefficient in the outer product
2321 VectorCoefficient * GetBCoef() const { return b; }
2322
2323 /// Evaluate the matrix coefficient at @a ip.
2325 const IntegrationPoint &ip) override;
2326};
2327
2328/** @brief Matrix coefficient defined as -a k x k x, for a vector k and scalar a
2329
2330 This coefficient returns $a * (|k|^2 I - k \otimes k)$, where I is
2331 the identity matrix and $\otimes$ indicates the outer product. This
2332 can be evaluated for vectors of any dimension but in three
2333 dimensions it corresponds to computing the cross product with k twice.
2334*/
2336{
2337private:
2338 real_t aConst;
2339 Coefficient * a;
2341
2342 mutable Vector vk;
2343
2344public:
2347
2348 /// Set the time for internally stored coefficients
2349 void SetTime(real_t t) override;
2350
2351 /// Reset the scalar factor as a constant
2352 void SetAConst(real_t A) { a = NULL; aConst = A; }
2353 /// Return the scalar factor
2354 real_t GetAConst() const { return aConst; }
2355
2356 /// Reset the scalar factor
2357 void SetACoef(Coefficient &A) { a = &A; }
2358 /// Return the scalar factor
2359 Coefficient * GetACoef() const { return a; }
2360
2361 /// Reset the vector factor
2362 void SetKCoef(VectorCoefficient &K) { k = &K; }
2363 /// Return the vector factor
2364 VectorCoefficient * GetKCoef() const { return k; }
2365
2366 /// Evaluate the matrix coefficient at @a ip.
2368 const IntegrationPoint &ip) override;
2369};
2370///@}
2371
2372/** @brief Vector quadrature function coefficient which requires that the
2373 quadrature rules used for this vector coefficient be the same as those that
2374 live within the supplied QuadratureFunction. */
2376{
2377private:
2378 const QuadratureFunction &QuadF; //do not own
2379 int index;
2380
2381public:
2382 /// Constructor with a quadrature function as input
2384
2385 /** Set the starting index within the QuadFunc that'll be used to project
2386 outwards as well as the corresponding length. The projected length should
2387 have the bounds of 1 <= length <= (length QuadFunc - index). */
2388 void SetComponent(int index_, int length_);
2389
2390 const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2391
2393 void Eval(Vector &V, ElementTransformation &T,
2394 const IntegrationPoint &ip) override;
2395
2396 void Project(QuadratureFunction &qf) override;
2397
2399};
2400
2401/** @brief Quadrature function coefficient which requires that the quadrature
2402 rules used for this coefficient be the same as those that live within the
2403 supplied QuadratureFunction. */
2405{
2406private:
2407 const QuadratureFunction &QuadF;
2408
2409public:
2410 /// Constructor with a quadrature function as input
2412
2413 const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2414
2415 real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override;
2416
2417 void Project(QuadratureFunction &qf) override;
2418
2420};
2421
2422/// Flags that determine what storage optimizations to use in CoefficientVector
2423enum class CoefficientStorage : int
2424{
2425 FULL = 0, ///< Store the coefficient as a full QuadratureFunction.
2426 CONSTANTS = 1 << 0, ///< Store constants using only @a vdim entries.
2427 SYMMETRIC = 1 << 1, ///< Store the triangular part of symmetric matrices.
2428 COMPRESSED = CONSTANTS | SYMMETRIC ///< Enable all above compressions.
2429};
2430
2435
2437{
2438 return int(a) & int(b);
2439}
2440
2441
2442/// @brief Class to represent a coefficient evaluated at quadrature points.
2443///
2444/// In the general case, a CoefficientVector is the same as a QuadratureFunction
2445/// with a coefficient projected onto it.
2446///
2447/// This class allows for some "compression" of the coefficient data, according
2448/// to the storage flags given by CoefficientStorage. For example, constant
2449/// coefficients can be stored using only @a vdim values, and symmetric matrices
2450/// can be stored using e.g. the upper triangular part of the matrix.
2452{
2453protected:
2454 CoefficientStorage storage; ///< Storage optimizations (see CoefficientStorage).
2455 int vdim; ///< Number of values per quadrature point.
2456 QuadratureSpaceBase &qs; ///< Associated QuadratureSpaceBase.
2457 QuadratureFunction *qf; ///< Internal QuadratureFunction (owned, may be NULL).
2458public:
2459 /// Create an empty CoefficientVector.
2462
2463 /// @brief Create a CoefficientVector from the given Coefficient and
2464 /// QuadratureSpaceBase.
2465 ///
2466 /// If @a coeff is NULL, it will be interpreted as a constant with value one.
2467 /// @sa CoefficientStorage for a description of @a storage_.
2470
2471 /// @brief Create a CoefficientVector from the given Coefficient and
2472 /// QuadratureSpaceBase.
2473 ///
2474 /// @sa CoefficientStorage for a description of @a storage_.
2477
2478 /// @brief Create a CoefficientVector from the given VectorCoefficient and
2479 /// QuadratureSpaceBase.
2480 ///
2481 /// @sa CoefficientStorage for a description of @a storage_.
2484
2485 /// @brief Create a CoefficientVector from the given MatrixCoefficient and
2486 /// QuadratureSpaceBase.
2487 ///
2488 /// @sa CoefficientStorage for a description of @a storage_.
2491
2492 /// @brief Evaluate the given Coefficient at the quadrature points defined by
2493 /// @ref qs.
2494 void Project(Coefficient &coeff);
2495
2496 /// @brief Evaluate the given VectorCoefficient at the quadrature points
2497 /// defined by @ref qs.
2498 ///
2499 /// @sa CoefficientVector for a description of the @a compress argument.
2500 void Project(VectorCoefficient &coeff);
2501
2502 /// @brief Evaluate the given MatrixCoefficient at the quadrature points
2503 /// defined by @ref qs.
2504 ///
2505 /// @sa CoefficientVector for a description of the @a compress argument.
2506 void Project(MatrixCoefficient &coeff, bool transpose=false);
2507
2508 /// @brief Project the transpose of @a coeff.
2509 ///
2510 /// @sa Project(MatrixCoefficient&, QuadratureSpace&, bool, bool)
2512
2513 /// Make this vector a reference to the given QuadratureFunction.
2514 void MakeRef(const QuadratureFunction &qf_);
2515
2516 /// Set this vector to the given constant.
2517 void SetConstant(real_t constant);
2518
2519 /// Set this vector to the given constant vector.
2520 void SetConstant(const Vector &constant);
2521
2522 /// Set this vector to the given constant matrix.
2523 void SetConstant(const DenseMatrix &constant, bool transpose=false);
2524
2525 /// Set this vector to the given constant symmetric matrix.
2526 void SetConstant(const DenseSymmetricMatrix &constant);
2527
2528 /// Return the number of values per quadrature point.
2529 int GetVDim() const;
2530
2532};
2533
2534/** @brief Compute the Lp norm of a function f.
2535 $ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} $ */
2537 const IntegrationRule *irs[]);
2538
2539/** @brief Compute the Lp norm of a vector function f = {f_i}_i=1...N.
2540 $ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} $ */
2542 const IntegrationRule *irs[]);
2543
2544#ifdef MFEM_USE_MPI
2545/** @brief Compute the global Lp norm of a function f.
2546 $ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} $ */
2548 const IntegrationRule *irs[]);
2549
2550/** @brief Compute the global Lp norm of a vector function f = {f_i}_i=1...N.
2551 $ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} $ */
2553 const IntegrationRule *irs[]);
2554#endif
2555
2556}
2557
2558#endif
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:1042
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...
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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:32
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.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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:65
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.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf with the piecewise constant values.
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 Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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:32
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.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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 Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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:234
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
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:46
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)