12 #ifndef MFEM_TEMPLATE_BILINEAR_FORM 13 #define MFEM_TEMPLATE_BILINEAR_FORM 15 #include "../config/tconfig.hpp" 16 #include "../linalg/simd.hpp" 17 #include "../linalg/ttensor.hpp" 41 template <
typename meshType,
typename solFESpace,
43 typename solVecLayout_t = ScalarLayout,
44 typename complex_t = double,
typename real_t = double,
45 typename impl_traits_t = AutoSIMDTraits<complex_t,real_t> >
62 static const int sdim = meshType::space_dim;
63 static const int dofs = solFE_type::dofs;
64 static const int vdim = solVecLayout_t::vec_dim;
65 static const int qpts = IR::qpts;
66 static const int AB = impl_traits_t::align_bytes;
67 static const int SS = impl_traits_t::simd_size;
68 static const int BE = impl_traits_t::batch_size;
72 typedef typename impl_traits_t::vreal_t
vreal_t;
78 typedef typename integ_t::coefficient_type
coeff_t;
80 typedef typename integ_t::template kernel<sdim,dim,vcomplex_t>::type
kernel_t;
82 typedef typename kernel_t::template p_asm_data<qpts>::type
p_assembled_t;
84 typedef typename kernel_t::template f_asm_data<qpts>::type
f_assembled_t;
87 typedef typename kernel_t::template
95 Trans_t::template Get<coeff_t,kernel_t>::EvalOps;
96 typedef typename Trans_t::template Result<EvalOps,impl_traits_t>
Type;
134 sol_fe(*sol_fes.FEColl()),
186 const int NE =
mesh.GetNE();
187 for (
int el = 0; el < NE; el +=
TE)
191 solFEval.Eval(el, R);
200 solFEval.Eval(el, R);
203 typename coeff_eval_t::result_t res;
206 for (
int k = 0; k <
BE; k++)
208 kernel_t::Action(k, F, wQ, res, R);
211 solFEval.template Assemble<true>(R);
221 const int NE =
mesh.GetNE();
224 const int size = ((NE+
TE-1)/
TE)*
BE;
227 for (
int el = 0; el < NE; el +=
TE)
232 typename coeff_eval_t::result_t res;
235 for (
int k = 0; k <
BE; k++)
242 inline MFEM_ALWAYS_INLINE
246 solFEval.
Eval(el, R);
248 for (
int k = 0; k <
BE; k++)
253 solFEval.template Assemble<true>(R);
264 const int NE =
mesh.GetNE();
265 for (
int el = 0; el < NE; el +=
TE)
271 #ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE 278 solFESpace solFESLoc(this->solFES);
282 const int NE =
mesh.GetNE();
283 for (
int el = 0; el < NE; el +=
TE)
285 solFESLoc.SetElement(el);
287 solFESLoc.VectorExtract(solVecLayoutLoc, x, xy_dof.
layout, xy_dof);
288 solFESLoc.VectorAssemble(xy_dof.
layout, xy_dof, solVecLayoutLoc, y);
295 typedef typename meshType::FESpace_type meshFESpace;
296 meshFESpace meshFES(
mesh.t_fes);
299 const int NE =
mesh.GetNE();
302 const int NVE = (NE+
TE-1)/
TE;
306 for (
int el = 0; el < NE; el +=
TE)
308 meshFES.SetElement(el);
309 meshFES.VectorExtract(
mesh.node_layout,
mesh.Nodes,
310 lnodes_t::layout, vsNodes);
311 vsNodes += lnodes_t::size;
322 const int NE =
mesh.GetNE();
325 const int size = ((NE+
TE-1)/
TE)*
BE;
329 for (
int el = 0; el < NE; el +=
TE)
332 T.EvalSerialized(el, vsNodes, F);
334 typename coeff_eval_t::result_t res;
337 for (
int k = 0; k <
BE; k++)
349 solVecLayout_t solVecLayoutLoc(this->solVecLayout);
350 solFESpace solFESLoc(this->solFES);
352 const int NE =
mesh.GetNE();
355 const int NVE = (NE+
TE-1)/
TE;
359 for (
int el = 0; el < NE; el +=
TE)
361 solFESLoc.SetElement(el);
362 solFESLoc.VectorExtract(solVecLayoutLoc, x, vdof_data_t::layout, vsx);
363 vsx += vdof_data_t::size;
373 const int NE =
mesh.GetNE();
377 for (
int el = 0; el < NE; el +=
TE)
382 for (
int k = 0; k <
BE; k++)
387 solFEval.template AssembleSerialized<false>(R, vsy);
393 #endif // MFEM_TEMPLATE_ENABLE_SERIALIZE 400 solFESpace solFESLoc(this->solFES);
402 solVecLayout_t solVecLayoutLoc(this->solVecLayout);
405 const int NE =
mesh.GetNE();
406 for (
int el = 0; el < NE; el +=
TE)
413 typename coeff_eval_t::result_t res;
416 for (
int k = 0; k <
BE; k++)
418 kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
424 for (
int k = 0; k <
BE; k++)
426 const int el_k = el+
SS*k;
427 if (el_k >= NE) {
break; }
430 S_spec::ElementMatrix::Compute(
431 asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.
layout, M_loc,
434 solFESLoc.SetElement(el_k);
435 for (
int bi = 0; bi <
vdim; bi++)
437 solFESLoc.AssembleBlock(bi, bi, solVecLayoutLoc, M_loc, M);
451 const int NE =
mesh.GetNE();
452 for (
int el = 0; el < NE; el +=
TE)
459 typename coeff_eval_t::result_t res;
462 for (
int k = 0; k <
BE; k++)
464 kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
471 for (
int k = 0; k <
BE; k++)
473 const int el_k = el+
SS*k;
474 if (el_k >= NE) {
break; }
477 S_spec::ElementMatrix::Compute(
478 asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.
layout, M_loc,
481 for (
int s = 0;
s <
SS && el_k+
s < NE;
s++)
483 complex_t *M_data = M.
GetData(el_k+
s);
484 for (
int j = 0; j <
dofs; j++)
486 for (
int i = 0; i <
dofs; i++)
488 M_data[j+
dofs*i] = M_loc(i,j)[
s];
506 const int *dof_map_ = dof_map->
GetData();
509 const int NE =
mesh.GetNE();
510 for (
int el = 0; el < NE; el +=
TE)
517 typename coeff_eval_t::result_t res;
520 for (
int k = 0; k <
BE; k++)
522 kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
528 for (
int k = 0; k <
BE; k++)
530 const int el_k = el+
SS*k;
531 if (el_k >= NE) {
break; }
534 S_spec::ElementMatrix::Compute(
535 asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.
layout, M_loc,
540 for (
int s = 0;
s <
SS && el_k+
s < NE;
s++)
542 for (
int i = 0; i <
dofs; i++)
544 for (
int j = 0; j <
dofs; j++)
546 M_loc_perm(dof_map_[i],dof_map_[j]) = M_loc(i,j)[
s];
549 for (
int bi = 1; bi <
vdim; bi++)
554 a.AssembleElementMatrix(el_k+
s, M_loc_perm, vdofs);
562 a.AssembleElementMatrix(el_k, DM, vdofs);
566 for (
int bi = 0; bi <
vdim; bi++)
570 a.AssembleElementMatrix(el_k, M_loc_perm, vdofs);
575 for (
int s = 0;
s <
SS && el_k+
s < NE;
s++)
577 for (
int i = 0; i <
dofs; i++)
579 for (
int j = 0; j <
dofs; j++)
581 M_loc_perm(i,j) = M_loc(i,j)[
s];
584 for (
int bi = 1; bi <
vdim; bi++)
589 a.AssembleElementMatrix(el_k+
s, M_loc_perm, vdofs);
603 solVecLayout_t solVecLayoutLoc(this->solVecLayout);
604 const int NE =
mesh.GetNE();
605 for (
int el = 0; el < NE; el++)
610 solFES.VectorExtract(solVecLayoutLoc, x, x_dof.
layout, x_dof);
613 x_dof.
layout.merge_23(), x_dof,
614 y_dof.
layout.merge_23(), y_dof);
615 solFES.VectorAssemble(y_dof.
layout, y_dof, solVecLayoutLoc, y);
623 #endif // MFEM_TEMPLATE_BILINEAR_FORM
Host memory; aligned at 64 bytes.
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
void Delete()
Delete the owned pointers and reset the Memory object.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Operator application: y+=A(x) (default) or y+=a*A(x).
static const layout_type layout
MFEM_ALWAYS_INLINE void Eval(DataType &F)
T * GetData()
Returns the data.
Data type dense matrix using column-major storage.
MFEM_ALWAYS_INLINE void EvalSerialized(const typename DataType::vcomplex_t *loc_dofs, DataType &F)
static const layout_type layout
Host memory; aligned at 32 bytes.
bool Empty() const
Return true if the Memory object is empty, see Reset().
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
complex_t - dof/qpt data type, real_t - ShapeEvaluator (FE basis) data type
void MakeDataOwner() const
Set the Vector data (host pointer) ownership flag.
Host memory; using new[] and delete[].
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
data_t data[aligned_size >0?aligned_size:1]
Rank 3 tensor (array of matrices)
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...