12 #ifndef MFEM_LIBCEED_COEFF
13 #define MFEM_LIBCEED_COEFF
17 #include "../../../general/forall.hpp"
18 #include "../../../config/config.hpp"
19 #include "../../../linalg/vector.hpp"
20 #include "../../../linalg/dtensor.hpp"
21 #include "../../../mesh/mesh.hpp"
22 #include "../../gridfunc.hpp"
23 #include "../../qfunction.hpp"
31 class IntegrationRule;
33 class VectorCoefficient;
53 virtual bool IsConstant()
const override {
return false; }
90 template <
typename Context>
99 coeff_ptr = ceedCoeff;
102 dynamic_cast<ConstantCoefficient*>(Q))
105 ctx.coeff = const_coeff->constant;
106 coeff_ptr = ceedCoeff;
109 dynamic_cast<GridFunctionCoefficient*>(Q))
113 coeff_ptr = ceedCoeff;
116 dynamic_cast<QuadratureFunctionCoefficient*>(Q))
119 const int ne = mesh.
GetNE();
122 MFEM_VERIFY(qFun.
Size() == nq * ne,
123 "Incompatible QuadratureFunction dimension \n");
126 "IntegrationRule used within integrator and in"
127 " QuadratureFunction appear to be different");
129 ceedCoeff->
coeff.
MakeRef(const_cast<mfem::QuadratureFunction &>(qFun),0);
131 coeff_ptr = ceedCoeff;
136 const int ne = mesh.
GetNE();
140 for (
int e = 0; e < ne; ++e)
143 for (
int q = 0; q < nq; ++q)
149 coeff_ptr = ceedCoeff;
164 template <
typename Context>
170 dynamic_cast<VectorConstantCoefficient*>(VQ))
172 const int vdim = const_coeff->GetVDim();
175 for (
int i = 0; i < vdim; i++)
177 ctx.coeff[i] = val[i];
179 coeff_ptr = ceedCoeff;
182 dynamic_cast<VectorGridFunctionCoefficient*>(VQ))
186 coeff_ptr = ceedCoeff;
189 dynamic_cast<VectorQuadratureFunctionCoefficient*>(VQ))
193 const int ne = mesh.
GetNE();
196 MFEM_VERIFY(qFun.
Size() == dim * nq * ne,
197 "Incompatible QuadratureFunction dimension \n");
200 "IntegrationRule used within integrator and in"
201 " QuadratureFunction appear to be different");
203 ceedCoeff->
coeff.
MakeRef(const_cast<mfem::QuadratureFunction &>(qFun),0);
205 coeff_ptr = ceedCoeff;
211 const int ne = mesh.
GetNE();
216 for (
int e = 0; e < ne; ++e)
219 VQ->
Eval(Q_ir, T, ir);
220 for (
int q = 0; q < nq; ++q)
222 for (
int i = 0; i <
dim; ++i)
224 C(i,q,e) = Q_ir(i,q);
229 coeff_ptr = ceedCoeff;
246 template <
typename Context>
257 coeff_ptr = ceedCoeff;
260 dynamic_cast<ConstantCoefficient*>(Q))
263 ctx.coeff = const_coeff->constant;
264 coeff_ptr = ceedCoeff;
267 dynamic_cast<GridFunctionCoefficient*>(Q))
271 coeff_ptr = ceedCoeff;
274 dynamic_cast<QuadratureFunctionCoefficient*>(Q))
277 const int ne = mesh.
GetNE();
280 MFEM_VERIFY(qFun.
Size() == nq * ne,
281 "Incompatible QuadratureFunction dimension \n");
284 "IntegrationRule used within integrator and in"
285 " QuadratureFunction appear to be different");
287 Memory<int> m_indices((
int*)indices, nelem,
false);
289 auto d_indices =
Read(m_indices, nelem);
291 MFEM_FORALL(i, nelem * nq,
294 const int sub_e = i/nq;
295 const int e = d_indices[sub_e];
296 out(q, sub_e) = in(q, e);
300 coeff_ptr = ceedCoeff;
308 for (
int i = 0; i < nelem; ++i)
310 const int e = indices[i];
312 for (
int q = 0; q < nq; ++q)
318 coeff_ptr = ceedCoeff;
336 template <
typename Context>
344 dynamic_cast<VectorConstantCoefficient*>(VQ))
346 const int vdim = const_coeff->GetVDim();
349 for (
int i = 0; i < vdim; i++)
351 ctx.coeff[i] = val[i];
353 coeff_ptr = ceedCoeff;
356 dynamic_cast<VectorGridFunctionCoefficient*>(VQ))
360 coeff_ptr = ceedCoeff;
363 dynamic_cast<VectorQuadratureFunctionCoefficient*>(VQ))
367 const int ne = mesh.
GetNE();
370 MFEM_VERIFY(qFun.
Size() == dim * nq * ne,
371 "Incompatible QuadratureFunction dimension \n");
374 "IntegrationRule used within integrator and in"
375 " QuadratureFunction appear to be different");
377 Memory<int> m_indices((
int*)indices, nelem,
false);
379 auto d_indices =
Read(m_indices, nelem);
381 MFEM_FORALL(i, nelem * nq,
384 const int sub_e = i/nq;
385 const int e = d_indices[sub_e];
386 for (
int d = 0; d <
dim; d++)
388 out(d, q, sub_e) = in(d, q, e);
393 coeff_ptr = ceedCoeff;
403 for (
int i = 0; i < nelem; ++i)
405 const int e = indices[i];
407 VQ->
Eval(Q_ir, T, ir);
408 for (
int q = 0; q < nq; ++q)
410 for (
int d = 0; d <
dim; ++d)
412 C(d, q, i) = Q_ir(d, q);
417 coeff_ptr = ceedCoeff;
421 template <
typename Coeff,
typename Context>
442 #endif // MFEM_LIBCEED_COEFF
int GetNPoints() const
Returns the number of the points in the integration rule.
Class for an integration rule - an Array of IntegrationPoint.
QuadCoefficient(int ncomp_)
Class for grid function - Vector with associated FE space.
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
void InitVector(const mfem::Vector &v, CeedVector &cv)
Initialize a CeedVector from an mfem::Vector.
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 ...
void SetSize(int s)
Resize the vector to size s.
CeedElemRestriction restr
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
Vector coefficient that is constant in space and time.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Data type dense matrix using column-major storage.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
const mfem::GridFunction & gf
GridCoefficient(const mfem::GridFunction &gf_)
virtual bool IsConstant() const
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
void InitCoefficient(mfem::Coefficient *Q, mfem::Mesh &mesh, const mfem::IntegrationRule &ir, Coefficient *&coeff_ptr, Context &ctx)
Initializes an mfem::ceed::Coefficient coeff_ptr from an mfem::Coefficient Q, an mfem::Mesh mesh...
struct s_NavierContext ctx
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
VariableCoefficient(int ncomp_, CeedEvalMode emode_)
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
CeedElemRestriction restr
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
void DeleteDevice(bool copy_to_host=true)
Delete the device pointer, if owned. If copy_to_host is true and the data is valid only on device...
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
virtual bool IsConstant() const override
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void InitCoefficientWithIndices(mfem::Coefficient *Q, mfem::Mesh &mesh, const mfem::IntegrationRule &ir, int nelem, const int *indices, Coefficient *&coeff_ptr, Context &ctx)
Initializes an mfem::ceed::Coefficient coeff_ptr from an mfem::Coefficient Q, an mfem::Mesh mesh...
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Vector coefficient defined by a vector GridFunction.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Represents values or vectors of values at quadrature points on a mesh.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.