MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
restriction.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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_LIBCEED_RESTR
13#define MFEM_LIBCEED_RESTR
14
15#include "ceed.hpp"
16
17namespace mfem
18{
19
20namespace ceed
21{
22
23#ifdef MFEM_USE_CEED
24/** @brief Initialize a CeedElemRestriction for non-mixed meshes.
25
26 @param[in] fes Input finite element space.
27 @param[in] ceed Input Ceed object.
28 @param[out] restr The address of the initialized CeedElemRestriction object.
29*/
30void InitRestriction(const FiniteElementSpace &fes,
31 Ceed ceed,
32 CeedElemRestriction *restr);
33
34/** @brief Initialize a CeedElemRestriction for mixed meshes.
35
36 @param[in] fes The finite element space.
37 @param[in] ceed The Ceed object.
38 @param[in] nelem The number of elements.
39 @param[in] indices The indices of the elements of same type in the
40 `FiniteElementSpace`.
41 @param[out] restr The `CeedElemRestriction` to initialize. */
42void InitRestrictionWithIndices(const FiniteElementSpace &fes,
43 int nelem,
44 const int* indices,
45 Ceed ceed,
46 CeedElemRestriction *restr);
47
48/** @brief Initialize a strided CeedElemRestriction
49
50 @param[in] fes Input finite element space.
51 @param[in] nelem is the number of elements.
52 @param[in] nqpts is the total number of quadrature points.
53 @param[in] qdatasize is the number of data per quadrature point.
54 @param[in] strides Array for strides between [nodes, components, elements].
55 Data for node i, component j, element k can be found in the L-vector at
56 index i*strides[0] + j*strides[1] + k*strides[2]. CEED_STRIDES_BACKEND may
57 be used with vectors created by a Ceed backend.
58 @param[out] restr The `CeedElemRestriction` to initialize. */
60 CeedInt nelem, CeedInt nqpts, CeedInt qdatasize,
61 const CeedInt *strides,
62 CeedElemRestriction *restr);
63
64/** @brief Initialize a CeedElemRestriction for a mfem::Coefficient on a mixed
65 mesh.
66
67 @param[in] fes The finite element space.
68 @param[in] nelem is the number of elements.
69 @param[in] indices The indices of the elements of same type in the
70 `FiniteElementSpace`.
71 @param[in] nquads is the total number of quadrature points
72 @param[in] ncomp is the number of data per quadrature point
73 @param[in] ceed The Ceed object.
74 @param[out] restr The `CeedElemRestriction` to initialize. */
75void InitCoeffRestrictionWithIndices(const FiniteElementSpace &fes,
76 int nelem,
77 const int* indices,
78 int nquads,
79 int ncomp,
80 Ceed ceed,
81 CeedElemRestriction *restr);
82
83#endif
84
85} // namespace ceed
86
87} // namespace mfem
88
89#endif // MFEM_LIBCEED_RESTR
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
void InitCoeffRestrictionWithIndices(const FiniteElementSpace &fes, int nelem, const int *indices, int nquads, int ncomp, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for a mfem::Coefficient on a mixed mesh.
void InitStridedRestriction(const mfem::FiniteElementSpace &fes, CeedInt nelem, CeedInt nqpts, CeedInt qdatasize, const CeedInt *strides, CeedElemRestriction *restr)
Initialize a strided CeedElemRestriction.
void InitRestriction(const FiniteElementSpace &fes, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for non-mixed meshes.
void InitRestrictionWithIndices(const FiniteElementSpace &fes, int nelem, const int *indices, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for mixed meshes.