12 #ifndef MFEM_TEMPLATE_BILINEAR_FORM
13 #define MFEM_TEMPLATE_BILINEAR_FORM
15 #include "../config/tconfig.hpp"
16 #include "../linalg/ttensor.hpp"
30 template <
typename meshType,
typename solFESpace,
31 typename IR,
typename IntegratorType,
32 typename solVecLayout_t = ScalarLayout,
33 typename complex_t = double,
typename real_t =
double>
47 static const int sdim = meshType::space_dim;
48 static const int dofs = solFE_type::dofs;
49 static const int vdim = solVecLayout_t::vec_dim;
50 static const int qpts = IR::qpts;
53 typedef typename integ_t::coefficient_type
coeff_t;
54 typedef typename integ_t::template kernel<sdim,dim,complex_t>::type
kernel_t;
55 typedef typename kernel_t::template p_asm_data<qpts>::type
p_assembled_t;
56 typedef typename kernel_t::template f_asm_data<qpts>::type
f_assembled_t;
62 Trans_t::template Get<coeff_t,kernel_t>::EvalOps;
63 typedef typename Trans_t::template Result<EvalOps,NE>
Type;
94 mesh(*sol_fes.GetMesh()),
114 const int num_elem = 1;
115 MultAssembled<num_elem>(x, y);
129 typedef typename kernel_t::template
130 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
141 const int NE =
mesh.GetNE();
142 for (
int el = 0; el < NE; el++)
146 solFEval.Eval(el, R);
155 solFEval.Eval(el, R);
158 typename coeff_eval_t::result_t res;
161 kernel_t::Action(0, F, wQ, res, R);
163 solFEval.template Assemble<true>(R);
171 typedef typename kernel_t::template
172 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
177 const int NE =
mesh.GetNE();
182 for (
int el = 0; el < NE; el++)
187 typename coeff_eval_t::result_t res;
190 for (
int k = 0; k < BE; k++)
197 template <
int num_elem>
198 inline MFEM_ALWAYS_INLINE
202 solFEval.
Eval(el, R);
204 for (
int k = 0; k < num_elem; k++)
209 solFEval.template Assemble<true>(R);
213 template <
int num_elem>
221 const int NE =
mesh.GetNE();
222 const int bNE = NE-NE%num_elem;
223 for (
int el = 0; el < bNE; el += num_elem)
225 ElementAddMultAssembled<num_elem>(el, solFEval);
227 for (
int el = bNE; el < NE; el++)
229 ElementAddMultAssembled<1>(el, solFEval);
233 #ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
244 const int NE =
mesh.GetNE();
245 for (
int el = 0; el < NE; el++)
247 solFES.SetElement(el);
249 solFES.VectorExtract(solVecLayout, x, xy_dof.
layout, xy_dof);
250 solFES.VectorAssemble(xy_dof.
layout, xy_dof, solVecLayout, y);
257 typedef typename meshType::FESpace_type meshFESpace;
258 meshFESpace meshFES(
mesh.t_fes);
261 const int NE =
mesh.GetNE();
262 sNodes.
SetSize(lnodes_t::size*NE);
264 for (
int el = 0; el < NE; el++)
266 meshFES.SetElement(el);
267 meshFES.VectorExtract(
mesh.node_layout,
mesh.Nodes,
268 lnodes_t::layout, lNodes);
269 lNodes += lnodes_t::size;
278 typedef typename kernel_t::template
279 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
284 const int NE =
mesh.GetNE();
289 for (
int el = 0; el < NE; el++)
294 typename coeff_eval_t::result_t res;
308 const int NE =
mesh.GetNE();
310 complex_t *loc_sx = sx.
GetData();
311 for (
int el = 0; el < NE; el++)
313 solFES.SetElement(el);
314 solFES.VectorExtract(solVecLayout, x, vdof_data_t::layout, loc_sx);
325 const int NE =
mesh.GetNE();
326 const complex_t *loc_sx = sx.
GetData();
327 complex_t *loc_sy = sy.
GetData();
328 for (
int el = 0; el < NE; el++)
335 solFEval.template AssembleSerialized<false>(R, loc_sy);
341 #endif // MFEM_TEMPLATE_ENABLE_SERIALIZE
348 typedef typename kernel_t::template
349 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
357 const int NE =
mesh.GetNE();
358 for (
int el = 0; el < NE; el++)
365 typename coeff_eval_t::result_t res;
368 kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
375 asm_qpt_data.layout, asm_qpt_data, M_loc.
layout, M_loc, solEval);
377 solFES.SetElement(el);
378 for (
int bi = 0; bi <
vdim; bi++)
380 solFES.AssembleBlock(bi, bi, solVecLayout, M_loc, M);
390 typedef typename kernel_t::template
391 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
397 const int NE =
mesh.GetNE();
398 for (
int el = 0; el < NE; el++)
405 typename coeff_eval_t::result_t res;
408 kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
416 asm_qpt_data.layout, asm_qpt_data, M_loc.
layout, M_loc, solEval);
418 complex_t *M_data = M.
GetData(el);
419 M_loc.template AssignTo<AssignOp::Set>(M_data);
428 typedef typename kernel_t::template
429 CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
437 const int *dof_map_ = dof_map->
GetData();
440 const int NE =
mesh.GetNE();
441 for (
int el = 0; el < NE; el++)
448 typename coeff_eval_t::result_t res;
451 kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
458 asm_qpt_data.layout, asm_qpt_data, M_loc.
layout, M_loc, solEval);
462 for (
int i = 0; i <
dofs; i++)
464 for (
int j = 0; j <
dofs; j++)
466 M_loc_perm(dof_map_[i],dof_map_[j]) = M_loc(i,j);
469 for (
int bi = 1; bi <
vdim; bi++)
471 M_loc_perm.
CopyMN(M_loc_perm, dofs, dofs, 0, 0,
485 for (
int bi = 0; bi <
vdim; bi++)
503 const int NE =
mesh.GetNE();
504 for (
int el = 0; el < NE; el++)
509 solFES.VectorExtract(solVecLayout, x, x_dof.
layout, x_dof);
512 x_dof.
layout.merge_23(), x_dof,
513 y_dof.
layout.merge_23(), y_dof);
514 solFES.VectorAssemble(y_dof.
layout, y_dof, solVecLayout, y);
521 #endif // MFEM_TEMPLATE_BILINEAR_FORM
void SetSize(int s)
Resize the vector if the new size is different.
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.
static const layout_type layout
MFEM_ALWAYS_INLINE void EvalSerialized(const complex_t *loc_dofs, DataType &F)
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)