22#if !defined(_WIN32) || !defined(_MSC_VER)
25#define stat(dir, buf) _stat(dir, buf)
26#define S_ISDIR(mode) (((mode) & _S_IFMT) == _S_IFDIR)
44 auto itb = mfem::internal::ceed_basis_map.begin();
45 while (itb != mfem::internal::ceed_basis_map.end())
47 if (std::get<0>(itb->first)==fes)
49 CeedBasisDestroy(&itb->second);
50 itb = mfem::internal::ceed_basis_map.erase(itb);
57 auto itr = mfem::internal::ceed_restr_map.begin();
58 while (itr != mfem::internal::ceed_restr_map.end())
60 if (std::get<0>(itr->first)==fes)
62 CeedElemRestrictionDestroy(&itr->second);
63 itr = mfem::internal::ceed_restr_map.erase(itr);
77 CeedVectorCreate(mfem::internal::ceed, v.
Size(), &cv);
80 CeedGetPreferredMemType(mfem::internal::ceed, &mem);
83 cv_ptr =
const_cast<CeedScalar*
>(v.
Read());
87 cv_ptr =
const_cast<CeedScalar*
>(v.
HostRead());
90 CeedVectorSetArray(cv, mem, CEED_USE_POINTER, cv_ptr);
95 Ceed ceed, CeedBasis *basis,
96 CeedElemRestriction *restr)
106 Ceed ceed, CeedBasis *basis,
107 CeedElemRestriction *restr)
117 Ceed ceed, CeedBasis *basis,
118 CeedElemRestriction *restr)
135 ierr = CeedOperatorGetCeed(oper, &ceed); PCeedChk(ierr);
139 ierr = CeedOperatorIsComposite(oper, &isComposite); PCeedChk(ierr);
140 CeedOperator *subops;
143#if CEED_VERSION_GE(0, 10, 2)
144 ierr = CeedCompositeOperatorGetSubList(oper, &subops); PCeedChk(ierr);
146 ierr = CeedOperatorGetSubList(oper, &subops); PCeedChk(ierr);
148 ierr = CeedOperatorGetQFunction(subops[0], &qf); PCeedChk(ierr);
152 ierr = CeedOperatorGetQFunction(oper, &qf); PCeedChk(ierr);
154 CeedInt numinputfields, numoutputfields;
155 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
156 CeedOperatorField *inputfields;
159 ierr = CeedOperatorGetFields(subops[0], &numinputfields, &inputfields,
160 &numoutputfields, NULL); PCeedChk(ierr);
164 ierr = CeedOperatorGetFields(oper, &numinputfields, &inputfields,
165 &numoutputfields, NULL); PCeedChk(ierr);
168 CeedVector if_vector;
170 int found_index = -1;
171 for (
int i = 0; i < numinputfields; ++i)
173 ierr = CeedOperatorFieldGetVector(inputfields[i], &if_vector); PCeedChk(ierr);
174 if (if_vector == CEED_VECTOR_ACTIVE)
178 return CeedError(ceed, 1,
"Multiple active vectors in CeedOperator!");
186 return CeedError(ceed, 1,
"No active vector in CeedOperator!");
188 *field = inputfields[found_index];
259 const char *install_dir = MFEM_INSTALL_DIR
"/include/mfem/fem/ceed";
260 const char *source_dir = MFEM_SOURCE_DIR
"/fem/ceed";
262 if (stat(install_dir, &m_stat) == 0 && S_ISDIR(m_stat.st_mode))
266 else if (stat(source_dir, &m_stat) == 0 && S_ISDIR(m_stat.st_mode))
272 MFEM_ABORT(
"Cannot find libCEED kernels in MFEM_INSTALL_DIR or "
static const IntegrationRule & GetRule(const FiniteElement &el, ElementTransformation &Trans)
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Abstract class for all finite elements.
Class for an integration rule - an Array of IntegrationPoint.
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
static const IntegrationRule & GetRule(const FiniteElement &fe, ElementTransformation &T)
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
int Size() const
Returns the size of the vector.
void trans(const Vector &u, Vector &x)
const IntegrationRule & GetRule< VectorDiffusionIntegrator >(const VectorDiffusionIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
const std::string & GetCeedPath()
Return the path to the libCEED q-function headers.
void InitVector(const mfem::Vector &v, CeedVector &cv)
Initialize a CeedVector from an mfem::Vector.
void RemoveBasisAndRestriction(const mfem::FiniteElementSpace *fes)
Remove from ceed_basis_map and ceed_restr_map the entries associated with the given fes.
const IntegrationRule & GetRule< VectorMassIntegrator >(const VectorMassIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
void InitBasisWithIndices(const FiniteElementSpace &fes, const IntegrationRule &ir, int nelem, const int *indices, Ceed ceed, CeedBasis *basis)
Initialize a CeedBasis for mixed meshes.
void InitRestriction(const FiniteElementSpace &fes, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for non-mixed meshes.
const IntegrationRule & GetRule< MassIntegrator >(const MassIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
void InitBasisAndRestrictionWithIndices(const FiniteElementSpace &fes, const IntegrationRule &irm, int nelem, const int *indices, Ceed ceed, CeedBasis *basis, CeedElemRestriction *restr)
void InitRestrictionWithIndices(const FiniteElementSpace &fes, int nelem, const int *indices, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for mixed meshes.
void InitBasisAndRestriction(const FiniteElementSpace &fes, const IntegrationRule &irm, Ceed ceed, CeedBasis *basis, CeedElemRestriction *restr)
Initialize a CeedBasis and a CeedElemRestriction based on an mfem::FiniteElementSpace fes,...
const IntegrationRule & GetRule< VectorConvectionNLFIntegrator >(const VectorConvectionNLFIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
int CeedOperatorGetActiveField(CeedOperator oper, CeedOperatorField *field)
const IntegrationRule & GetRule< ConvectionIntegrator >(const ConvectionIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
void InitBasis(const FiniteElementSpace &fes, const IntegrationRule &ir, Ceed ceed, CeedBasis *basis)
Initialize a CeedBasis for non-mixed meshes.
const IntegrationRule & GetRule< DiffusionIntegrator >(const DiffusionIntegrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &trans)
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
@ DEVICE_MASK
Biwise-OR of all device backends.
@ CEED_MASK
Bitwise-OR of all CEED backends.