30ConvectionIntegrator::Kernels::Kernels()
34 ConvectionIntegrator::AddSpecialization<2, 2, 2>();
35 ConvectionIntegrator::AddSpecialization<2, 3, 3>();
36 ConvectionIntegrator::AddSpecialization<2, 4, 4>();
37 ConvectionIntegrator::AddSpecialization<2, 5, 5>();
38 ConvectionIntegrator::AddSpecialization<2, 6, 6>();
40 ConvectionIntegrator::AddSpecialization<2, 2, 3>();
41 ConvectionIntegrator::AddSpecialization<2, 3, 4>();
42 ConvectionIntegrator::AddSpecialization<2, 4, 5>();
43 ConvectionIntegrator::AddSpecialization<2, 5, 6>();
44 ConvectionIntegrator::AddSpecialization<2, 6, 7>();
47 ConvectionIntegrator::AddSpecialization<3, 2, 2>();
48 ConvectionIntegrator::AddSpecialization<3, 3, 3>();
49 ConvectionIntegrator::AddSpecialization<3, 4, 4>();
50 ConvectionIntegrator::AddSpecialization<3, 5, 5>();
51 ConvectionIntegrator::AddSpecialization<3, 6, 6>();
53 ConvectionIntegrator::AddSpecialization<3, 2, 3>();
54 ConvectionIntegrator::AddSpecialization<3, 3, 4>();
55 ConvectionIntegrator::AddSpecialization<3, 4, 5>();
56 ConvectionIntegrator::AddSpecialization<3, 5, 6>();
57 ConvectionIntegrator::AddSpecialization<3, 6, 7>();
61static void PAConvectionSetup2D(
const int NQ,
63 const Array<real_t> &w,
69 constexpr int DIM = 2;
71 const bool const_v =
vel.Size() ==
DIM;
73 const auto W = w.Read();
75 const auto V = const_v ?
82 const int e = q_global / NQ;
83 const int q = q_global % NQ;
84 const real_t J11 = J(q,0,0,e);
85 const real_t J21 = J(q,1,0,e);
86 const real_t J12 = J(q,0,1,e);
87 const real_t J22 = J(q,1,1,e);
89 const real_t v0 = const_v ? V(0,0,0) : V(0,q,e);
90 const real_t v1 = const_v ? V(1,0,0) : V(1,q,e);
94 y(q,0,e) = wx * J22 - wy * J12;
95 y(q,1,e) = -wx * J21 + wy * J11;
100static void PAConvectionSetup3D(
const int NQ,
102 const Array<real_t> &w,
108 constexpr int DIM = 3;
110 const auto W =
Reshape(w.Read(), NQ);
112 const bool const_v =
vel.Size() ==
DIM;
113 const auto V = const_v ?
116 auto y =
Reshape(op.Write(), NQ,3,NE);
119 const int e = q_global / NQ;
120 const int q = q_global % NQ;
121 const real_t J11 = J(q,0,0,e);
122 const real_t J12 = J(q,0,1,e);
123 const real_t J13 = J(q,0,2,e);
124 const real_t J21 = J(q,1,0,e);
125 const real_t J22 = J(q,1,1,e);
126 const real_t J23 = J(q,1,2,e);
127 const real_t J31 = J(q,2,0,e);
128 const real_t J32 = J(q,2,1,e);
129 const real_t J33 = J(q,2,2,e);
131 const real_t v0 = const_v ? V(0,0,0) : V(0,q,e);
132 const real_t v1 = const_v ? V(1,0,0) : V(1,q,e);
133 const real_t v2 = const_v ? V(2,0,0) : V(2,q,e);
138 const real_t A11 = (J22 * J33) - (J23 * J32);
139 const real_t A12 = (J32 * J13) - (J12 * J33);
140 const real_t A13 = (J12 * J23) - (J22 * J13);
141 const real_t A21 = (J31 * J23) - (J21 * J33);
142 const real_t A22 = (J11 * J33) - (J13 * J31);
143 const real_t A23 = (J21 * J13) - (J11 * J23);
144 const real_t A31 = (J21 * J32) - (J31 * J22);
145 const real_t A32 = (J31 * J12) - (J11 * J32);
146 const real_t A33 = (J11 * J22) - (J12 * J21);
148 y(q,0,e) = wx * A11 + wy * A12 + wz * A13;
149 y(q,1,e) = wx * A21 + wy * A22 + wz * A23;
150 y(q,2,e) = wx * A31 + wy * A32 + wz * A33;
154static void PAConvectionSetup(
const int dim,
157 const Array<real_t> &W,
163 if (
dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAConvectionSetup"); }
166 PAConvectionSetup2D(NQ, NE, W, J, coeff,
alpha, op);
170 PAConvectionSetup3D(NQ, NE, W, J, coeff,
alpha, op);
174void ConvectionIntegrator::AssemblePA(
const FiniteElementSpace &fes)
176 const MemoryType mt = (pa_mt == MemoryType::DEFAULT) ?
177 Device::GetDeviceMemoryType() : pa_mt;
179 Mesh *mesh = fes.GetMesh();
180 const FiniteElement &el = *fes.GetTypicalFE();
181 ElementTransformation &Trans = *mesh->GetTypicalElementTransformation();
182 const IntegrationRule *ir = IntRule ? IntRule : &
GetRule(el, Trans);
186 const bool mixed = mesh->GetNumGeometries(mesh->Dimension()) > 1 ||
187 fes.IsVariableOrder();
190 ceedOp =
new ceed::MixedPAConvectionIntegrator(*
this, fes, Q,
alpha);
194 ceedOp =
new ceed::PAConvectionIntegrator(fes, *ir, Q,
alpha);
198 const int dims = el.GetDim();
199 const int symmDims = dims;
200 nq = ir->GetNPoints();
201 dim = mesh->Dimension();
203 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS, mt);
204 maps = &el.GetDofToQuad(*ir, DofToQuad::TENSOR);
207 pa_data.SetSize(symmDims * nq * ne, mt);
209 QuadratureSpace qs(*mesh, *ir);
210 CoefficientVector
vel(*Q, qs, CoefficientStorage::COMPRESSED);
212 PAConvectionSetup(
dim, nq, ne, ir->GetWeights(), geom->J,
216void ConvectionIntegrator::AssembleDiagonalPA(Vector &diag)
220 ceedOp->GetDiagonal(diag);
224 MFEM_ABORT(
"AssembleDiagonalPA not yet implemented for"
225 " ConvectionIntegrator.");
229inline ConvectionIntegrator::ApplyKernelType
230ConvectionIntegrator::ApplyPAKernels::Fallback(
int DIM,
int,
int)
234 return PAConvectionApply2D;
238 return PAConvectionApply3D;
246inline ConvectionIntegrator::ApplyKernelType
247ConvectionIntegrator::ApplyPATKernels::Fallback(
int DIM,
int,
int)
251 return PAConvectionApplyT2D;
255 return PAConvectionApplyT3D;
263void ConvectionIntegrator::AddMultPA(
const Vector &x, Vector &y)
const
267 ceedOp->AddMult(x, y);
271 ApplyPAKernels::Run(
dim, dofs1D, quad1D, ne, maps->B, maps->G, maps->Bt,
272 maps->Gt, pa_data, x, y, dofs1D, quad1D);
276void ConvectionIntegrator::AddMultTransposePA(
const Vector &x, Vector &y)
const
280 MFEM_ABORT(
"AddMultPA not yet implemented with libCEED for"
281 " ConvectionIntegrator.");
285 ApplyPATKernels::Run(
dim, dofs1D, quad1D, ne, maps->B, maps->G, maps->Bt,
286 maps->Gt, pa_data, x, y, dofs1D, quad1D);
ConvectionIntegrator(VectorCoefficient &q, real_t a=1.0)
const IntegrationRule & GetRule(const Integrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
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,...
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
SchrodingerBaseKernels< ParMesh, ParFiniteElementSpace, ParComplexGridFunction, ParGridFunction, ParBilinearForm, ParMixedBilinearForm, ParLinearForm > Kernels
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
MemoryType
Memory types supported by MFEM.
void forall(int N, lambda &&body)
void vel(const Vector &x, real_t t, Vector &u)