MFEM  v4.6.0
Finite element discretization library
tbilinearform.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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_TEMPLATE_BILINEAR_FORM
13 #define MFEM_TEMPLATE_BILINEAR_FORM
14 
15 #include "../config/tconfig.hpp"
16 #include "../linalg/simd.hpp"
17 #include "../linalg/ttensor.hpp"
18 #include "bilinearform.hpp"
19 #include "tevaluator.hpp"
20 #include "teltrans.hpp"
21 #include "tcoefficient.hpp"
22 #include "fespace.hpp"
23 
24 namespace mfem
25 {
26 
27 /** @brief Templated bilinear form class, cf. bilinearform.?pp
28 
29  @tparam meshType typically TMesh, which is templated on FE type
30  @tparam solFESpace eg. H1_FiniteElementSpace
31  @tparam IR integration rule, typically TIntegrationRule, which is further
32  templated on element geometry
33  @tparam IntegratorType typically a TIntegrator, which is templated on a
34  kernel, eg. TDiffusionKernel or TMassKernel. This
35  describes what actual problem you solve.
36  @tparam solVecLayout_t describes how degrees of freedom are laid out,
37  scalar or vector, column/row major, etc.
38  @tparam complex_t data type for solution dofs
39  @tparam real_t data type for mesh nodes, solution basis, and mesh basis
40 */
41 template <typename meshType, typename solFESpace,
42  typename IR, typename IntegratorType,
43  typename solVecLayout_t = ScalarLayout,
44  typename complex_t = double, typename real_t = double,
45  typename impl_traits_t = AutoSIMDTraits<complex_t,real_t> >
46 class TBilinearForm : public Operator
47 {
48 public:
49  typedef impl_traits_t impl_traits_type;
50 
51 protected:
52  typedef complex_t complex_type;
53  typedef real_t real_type;
54 
55  typedef typename meshType::FE_type meshFE_type;
57  typedef typename solFESpace::FE_type solFE_type;
59  typedef solVecLayout_t solVecLayout_type;
60 
61  static const int dim = meshType::dim;
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;
69  static const int TE = SS*BE;
70 
71  typedef typename impl_traits_t::vcomplex_t vcomplex_t;
72  typedef typename impl_traits_t::vreal_t vreal_t;
73 
74  /// @name IntegratorType defines several internal types
75  ///@{
77  /// coeff_t might be TConstantCoefficient or TFunctionCoefficient, for example
78  typedef typename integ_t::coefficient_type coeff_t;
79  /// kernel_t may be TDiffusionKernel or TMassKernel
80  typedef typename integ_t::template kernel<sdim,dim,vcomplex_t>::type kernel_t;
81  /// p_assembled_t is something like a TTensor or TMatrix for partial assembly
82  typedef typename kernel_t::template p_asm_data<qpts>::type p_assembled_t;
83  /// f_assembled_t is something like a TTensor or TMatrix for full assembly
84  typedef typename kernel_t::template f_asm_data<qpts>::type f_assembled_t;
85  ///@}
86 
87  typedef typename kernel_t::template
88  CoefficientEval<IR,coeff_t,impl_traits_t>::Type coeff_eval_t;
89 
90 
92  struct T_result
93  {
94  static const int EvalOps =
95  Trans_t::template Get<coeff_t,kernel_t>::EvalOps;
96  typedef typename Trans_t::template Result<EvalOps,impl_traits_t> Type;
97  };
98 
99  typedef FieldEvaluator<solFESpace,solVecLayout_t,IR,
100  complex_t,real_t> solFieldEval;
101 
102  /** @brief Contains matrix sizes, type of kernel (ElementMatrix is templated
103  on a kernel, e.g. ElementMatrix::Compute may be AssembleGradGrad()). */
104  struct S_spec
105  {
106  typedef typename solFieldEval::template Spec<kernel_t,impl_traits_t> Spec;
107  typedef typename Spec::DataType DataType;
108  typedef typename Spec::ElementMatrix ElementMatrix;
109  };
110 
111  // Data members
112 
113  meshType mesh;
115 
118  mutable solFESpace solFES;
119  solVecLayout_t solVecLayout;
120 
122 
124 
126 
128 
129 public:
130  TBilinearForm(const IntegratorType &integ, const FiniteElementSpace &sol_fes)
131  : Operator(sol_fes.GetNDofs()*vdim),
132  mesh(*sol_fes.GetMesh()),
133  meshEval(mesh.fe),
134  sol_fe(*sol_fes.FEColl()),
135  solEval(sol_fe),
136  solFES(sol_fe, sol_fes),
137  solVecLayout(sol_fes),
138  int_rule(),
139  coeff(integ.coeff),
140  assembled_data(),
141  in_fes(sol_fes)
142  {
144  AB == 32 ? MemoryType::HOST_32 :
146  }
147 
148  virtual ~TBilinearForm()
149  {
151  }
152 
153  /// Get the input finite element space prolongation matrix
154  virtual const Operator *GetProlongation() const
155  { return ((FiniteElementSpace &)in_fes).GetProlongationMatrix(); }
156  /// Get the input finite element space restriction matrix
157  virtual const Operator *GetRestriction() const
158  { return ((FiniteElementSpace &)in_fes).GetRestrictionMatrix(); }
159 
160  virtual void Mult(const Vector &x, Vector &y) const
161  {
162  if (!assembled_data.Empty())
163  {
164  MultAssembled(x, y);
165  }
166  else
167  {
168  MultUnassembled(x, y);
169  }
170  }
171 
172  // complex_t = double
173  void MultUnassembled(const Vector &x, Vector &y) const
174  {
175  y = 0.0;
176 
177  // For better performance, create stack copies of solFES, and solEval
178  // inside 'solFEval'. The element-transformation 'T' also copies the
179  // meshFES, meshEval, etc internally.
180  // Is performance actually better with this implementation?
181  Trans_t T(mesh, meshEval);
183  x.GetData(), y.GetData());
185 
186  const int NE = mesh.GetNE();
187  for (int el = 0; el < NE; el += TE)
188  {
189 #if 0
190  typename S_spec::DataType R;
191  solFEval.Eval(el, R);
192 
193  typename T_result::Type F;
194  T.Eval(el, F);
195 #else
196  typename T_result::Type F;
197  T.Eval(el, F);
198 
199  typename S_spec::DataType R;
200  solFEval.Eval(el, R);
201 #endif
202 
203  typename coeff_eval_t::result_t res;
204  wQ.Eval(F, res);
205 
206  for (int k = 0; k < BE; k++)
207  {
208  kernel_t::Action(k, F, wQ, res, R);
209  }
210 
211  solFEval.template Assemble<true>(R);
212  }
213  }
214 
215  /// Partial assembly of quadrature point data
216  void Assemble()
217  {
218  Trans_t T(mesh, meshEval);
220 
221  const int NE = mesh.GetNE();
222  if (assembled_data.Empty())
223  {
224  const int size = ((NE+TE-1)/TE)*BE;
226  }
227  for (int el = 0; el < NE; el += TE)
228  {
229  typename T_result::Type F;
230  T.Eval(el, F);
231 
232  typename coeff_eval_t::result_t res;
233  wQ.Eval(F, res);
234 
235  for (int k = 0; k < BE; k++)
236  {
237  kernel_t::Assemble(k, F, wQ, res, assembled_data[el/SS+k]);
238  }
239  }
240  }
241 
242  inline MFEM_ALWAYS_INLINE
243  void ElementAddMultAssembled(int el, solFieldEval &solFEval) const
244  {
245  typename S_spec::DataType R;
246  solFEval.Eval(el, R);
247 
248  for (int k = 0; k < BE; k++)
249  {
250  kernel_t::MultAssembled(k, assembled_data[el/SS+k], R);
251  }
252 
253  solFEval.template Assemble<true>(R);
254  }
255 
256  // complex_t = double
257  void MultAssembled(const Vector &x, Vector &y) const
258  {
259  y = 0.0;
260 
262  x.GetData(), y.GetData());
263 
264  const int NE = mesh.GetNE();
265  for (int el = 0; el < NE; el += TE)
266  {
267  ElementAddMultAssembled(el, solFEval);
268  }
269  }
270 
271 #ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
272  // complex_t = double
274  {
275  y = 0.0;
276 
277  solVecLayout_type solVecLayoutLoc(this->solVecLayout);
278  solFESpace solFESLoc(this->solFES);
279 
281 
282  const int NE = mesh.GetNE();
283  for (int el = 0; el < NE; el += TE)
284  {
285  solFESLoc.SetElement(el);
286 
287  solFESLoc.VectorExtract(solVecLayoutLoc, x, xy_dof.layout, xy_dof);
288  solFESLoc.VectorAssemble(xy_dof.layout, xy_dof, solVecLayoutLoc, y);
289  }
290  }
291 
292  // real_t = double
293  void SerializeNodes(Vector &sNodes) const
294  {
295  typedef typename meshType::FESpace_type meshFESpace;
296  meshFESpace meshFES(mesh.t_fes);
298 
299  const int NE = mesh.GetNE();
300  // TODO: How do we make sure that this array is aligned properly, AND the
301  // compiler knows that it is aligned? => ALIGN_32|ALIGN_64 when ready
302  const int NVE = (NE+TE-1)/TE;
303  vreal_t *vsNodes = new vreal_t[lnodes_t::size*NVE];
304  sNodes.NewDataAndSize(vsNodes[0].vec, (lnodes_t::size*SS)*NVE);
305  sNodes.MakeDataOwner();
306  for (int el = 0; el < NE; el += TE)
307  {
308  meshFES.SetElement(el);
309  meshFES.VectorExtract(mesh.node_layout, mesh.Nodes,
310  lnodes_t::layout, vsNodes);
311  vsNodes += lnodes_t::size;
312  }
313  }
314 
315  /// Partial assembly from "serialized" nodes
316  // real_t = double
318  {
319  Trans_t T(mesh, meshEval);
321 
322  const int NE = mesh.GetNE();
323  if (assembled_data.Empty())
324  {
325  const int size = ((NE+TE-1)/TE)*BE;
327  }
328  const vreal_t *vsNodes = (const vreal_t*)(sNodes.GetData());
329  for (int el = 0; el < NE; el += TE)
330  {
331  typename T_result::Type F;
332  T.EvalSerialized(el, vsNodes, F);
333 
334  typename coeff_eval_t::result_t res;
335  wQ.Eval(F, res);
336 
337  for (int k = 0; k < BE; k++)
338  {
339  kernel_t::Assemble(k, F, wQ, res, assembled_data[el/SS+k]);
340  }
341  }
342  }
343 
344  // complex_t = double
345  void Serialize(const Vector &x, Vector &sx) const
346  {
347  typedef TTensor3<dofs,vdim,BE,vcomplex_t> vdof_data_t;
348 
349  solVecLayout_t solVecLayoutLoc(this->solVecLayout);
350  solFESpace solFESLoc(this->solFES);
351 
352  const int NE = mesh.GetNE();
353  // TODO: How do we make sure that this array is aligned properly, AND
354  // the compiler knows that it is aligned? => ALIGN_32|ALIGN_64 when ready
355  const int NVE = (NE+TE-1)/TE;
356  vreal_t *vsx = new vreal_t[vdof_data_t::size*NVE];
357  sx.NewDataAndSize(vsx[0].vec, (vdof_data_t::size*SS)*NVE);
358  sx.MakeDataOwner();
359  for (int el = 0; el < NE; el += TE)
360  {
361  solFESLoc.SetElement(el);
362  solFESLoc.VectorExtract(solVecLayoutLoc, x, vdof_data_t::layout, vsx);
363  vsx += vdof_data_t::size;
364  }
365  }
366 
367  /// serialized vector sx --> serialized vector 'sy'
368  // complex_t = double
369  void MultAssembledSerialized(const Vector &sx, Vector &sy) const
370  {
371  solFieldEval solFEval(solFES, solEval, solVecLayout, NULL, NULL);
372 
373  const int NE = mesh.GetNE();
374  const vreal_t *vsx = (const vreal_t*)(sx.GetData());
375  vreal_t *vsy = (vreal_t*)(sy.GetData());
376 
377  for (int el = 0; el < NE; el += TE)
378  {
379  typename S_spec::DataType R;
380  solFEval.EvalSerialized(vsx, R);
381 
382  for (int k = 0; k < BE; k++)
383  {
384  kernel_t::MultAssembled(k, assembled_data[el/SS+k], R);
385  }
386 
387  solFEval.template AssembleSerialized<false>(R, vsy);
388 
389  vsx += vdim*dofs*BE;
390  vsy += vdim*dofs*BE;
391  }
392  }
393 #endif // MFEM_TEMPLATE_ENABLE_SERIALIZE
394 
395  /// Assemble the operator in a SparseMatrix.
396  // complex_t = double
398  {
399  Trans_t T(mesh, meshEval);
400  solFESpace solFESLoc(this->solFES);
401  solShapeEval solEvalLoc(this->solEval);
402  solVecLayout_t solVecLayoutLoc(this->solVecLayout);
404 
405  const int NE = mesh.GetNE();
406  for (int el = 0; el < NE; el += TE)
407  {
408  f_assembled_t asm_qpt_data[BE];
409  {
410  typename T_result::Type F;
411  T.Eval(el, F);
412 
413  typename coeff_eval_t::result_t res;
414  wQ.Eval(F, res);
415 
416  for (int k = 0; k < BE; k++)
417  {
418  kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
419  }
420  }
421 
422  // For now, when vdim > 1, assume block-diagonal matrix with the same
423  // diagonal block for all components.
424  for (int k = 0; k < BE; k++)
425  {
426  const int el_k = el+SS*k;
427  if (el_k >= NE) { break; }
428 
430  S_spec::ElementMatrix::Compute(
431  asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.layout, M_loc,
432  solEvalLoc);
433 
434  solFESLoc.SetElement(el_k);
435  for (int bi = 0; bi < vdim; bi++)
436  {
437  solFESLoc.AssembleBlock(bi, bi, solVecLayoutLoc, M_loc, M);
438  }
439  }
440  }
441  }
442 
443  /// Assemble element matrices and store them as a DenseTensor object.
444  // complex_t = double
446  {
447  Trans_t T(mesh, meshEval);
448  solShapeEval solEvalLoc(this->solEval);
450 
451  const int NE = mesh.GetNE();
452  for (int el = 0; el < NE; el += TE)
453  {
454  f_assembled_t asm_qpt_data[BE];
455  {
456  typename T_result::Type F;
457  T.Eval(el, F);
458 
459  typename coeff_eval_t::result_t res;
460  wQ.Eval(F, res);
461 
462  for (int k = 0; k < BE; k++)
463  {
464  kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
465  }
466  }
467 
468  // For now, when vdim > 1, assume block-diagonal matrix with the same
469  // diagonal block for all components.
470  // M is assumed to be (dof x dof x NE).
471  for (int k = 0; k < BE; k++)
472  {
473  const int el_k = el+SS*k;
474  if (el_k >= NE) { break; }
475 
477  S_spec::ElementMatrix::Compute(
478  asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.layout, M_loc,
479  solEvalLoc);
480 
481  for (int s = 0; s < SS && el_k+s < NE; s++)
482  {
483  complex_t *M_data = M.GetData(el_k+s);
484  for (int j = 0; j < dofs; j++)
485  {
486  for (int i = 0; i < dofs; i++)
487  {
488  M_data[j+dofs*i] = M_loc(i,j)[s];
489  }
490  }
491  }
492  }
493  }
494  }
495 
496  /// Assemble element matrices and add them to the bilinear form
497  // complex_t = double
499  {
500  Trans_t T(mesh, meshEval);
501  solShapeEval solEvalLoc(this->solEval);
503 
504  Array<int> vdofs;
505  const Array<int> *dof_map = sol_fe.GetDofMap();
506  const int *dof_map_ = dof_map->GetData();
507  DenseMatrix M_loc_perm(dofs*vdim,dofs*vdim); // initialized with zeros
508 
509  const int NE = mesh.GetNE();
510  for (int el = 0; el < NE; el += TE)
511  {
512  f_assembled_t asm_qpt_data[BE];
513  {
514  typename T_result::Type F;
515  T.Eval(el, F);
516 
517  typename coeff_eval_t::result_t res;
518  wQ.Eval(F, res);
519 
520  for (int k = 0; k < BE; k++)
521  {
522  kernel_t::Assemble(k, F, wQ, res, asm_qpt_data[k]);
523  }
524  }
525 
526  // For now, when vdim > 1, assume block-diagonal matrix with the same
527  // diagonal block for all components.
528  for (int k = 0; k < BE; k++)
529  {
530  const int el_k = el+SS*k;
531  if (el_k >= NE) { break; }
532 
534  S_spec::ElementMatrix::Compute(
535  asm_qpt_data[k].layout, asm_qpt_data[k], M_loc.layout, M_loc,
536  solEvalLoc);
537 
538  if (dof_map) // switch from tensor-product ordering
539  {
540  for (int s = 0; s < SS && el_k+s < NE; s++)
541  {
542  for (int i = 0; i < dofs; i++)
543  {
544  for (int j = 0; j < dofs; j++)
545  {
546  M_loc_perm(dof_map_[i],dof_map_[j]) = M_loc(i,j)[s];
547  }
548  }
549  for (int bi = 1; bi < vdim; bi++)
550  {
551  M_loc_perm.CopyMN(M_loc_perm, dofs, dofs, 0, 0,
552  bi*dofs, bi*dofs);
553  }
554  a.AssembleElementMatrix(el_k+s, M_loc_perm, vdofs);
555  }
556  }
557  else if (SS == 1)
558  {
559  DenseMatrix DM(M_loc.data[0].vec, dofs, dofs);
560  if (vdim == 1)
561  {
562  a.AssembleElementMatrix(el_k, DM, vdofs);
563  }
564  else
565  {
566  for (int bi = 0; bi < vdim; bi++)
567  {
568  M_loc_perm.CopyMN(DM, dofs, dofs, 0, 0, bi*dofs, bi*dofs);
569  }
570  a.AssembleElementMatrix(el_k, M_loc_perm, vdofs);
571  }
572  }
573  else
574  {
575  for (int s = 0; s < SS && el_k+s < NE; s++)
576  {
577  for (int i = 0; i < dofs; i++)
578  {
579  for (int j = 0; j < dofs; j++)
580  {
581  M_loc_perm(i,j) = M_loc(i,j)[s];
582  }
583  }
584  for (int bi = 1; bi < vdim; bi++)
585  {
586  M_loc_perm.CopyMN(M_loc_perm, dofs, dofs, 0, 0,
587  bi*dofs, bi*dofs);
588  }
589  a.AssembleElementMatrix(el_k+s, M_loc_perm, vdofs);
590  }
591  }
592  }
593  }
594  }
595 
596  /// Multiplication using assembled element matrices stored as a DenseTensor.
597  // complex_t = double
598  void AddMult(DenseTensor &M, const Vector &x, Vector &y) const
599  {
600  // For now, when vdim > 1, assume block-diagonal matrix with the same
601  // diagonal block for all components.
602  // M is assumed to be (dof x dof x NE).
603  solVecLayout_t solVecLayoutLoc(this->solVecLayout);
604  const int NE = mesh.GetNE();
605  for (int el = 0; el < NE; el++)
606  {
608 
609  solFES.SetElement(el);
610  solFES.VectorExtract(solVecLayoutLoc, x, x_dof.layout, x_dof);
611  Mult_AB<false>(TMatrix<dofs,dofs>::layout,
612  M(el).Data(),
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);
616  }
617  }
618  using Operator::AddMult;
619 };
620 
621 } // namespace mfem
622 
623 #endif // MFEM_TEMPLATE_BILINEAR_FORM
ShapeEvaluator< meshFE_type, IR, real_t > meshShapeEval
Templated bilinear form class, cf. bilinearform.?pp.
Host memory; aligned at 64 bytes.
Contains matrix sizes, type of kernel (ElementMatrix is templated on a kernel, e.g. ElementMatrix::Compute may be AssembleGradGrad()).
kernel_t::template p_asm_data< qpts >::type p_assembled_t
p_assembled_t is something like a TTensor or TMatrix for partial assembly
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:160
kernel_t::template f_asm_data< qpts >::type f_assembled_t
f_assembled_t is something like a TTensor or TMatrix for full assembly
static const int vdim
const FiniteElementSpace & in_fes
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).
Definition: operator.cpp:51
static const layout_type layout
Definition: ttensor.hpp:356
MFEM_ALWAYS_INLINE void Eval(DataType &F)
meshShapeEval meshEval
T * GetData()
Returns the data.
Definition: array.hpp:115
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Trans_t::template Result< EvalOps, impl_traits_t > Type
MFEM_ALWAYS_INLINE void EvalSerialized(const typename DataType::vcomplex_t *loc_dofs, DataType &F)
static const layout_type layout
Definition: ttensor.hpp:392
Element transformation class, templated on a mesh type and an integration rule. It is constructed fro...
Definition: teltrans.hpp:33
void SerializeNodes(Vector &sNodes) const
void AssembleMatrix(DenseTensor &M) const
Assemble element matrices and store them as a DenseTensor object.
Host memory; aligned at 32 bytes.
void AssembleBilinearForm(BilinearForm &a) const
Assemble element matrices and add them to the bilinear form.
integ_t::template kernel< sdim, dim, vcomplex_t >::type kernel_t
kernel_t may be TDiffusionKernel or TMassKernel
bool Empty() const
Return true if the Memory object is empty, see Reset().
Data type sparse matrix.
Definition: sparsemat.hpp:50
solVecLayout_t solVecLayout_type
solFieldEval::template Spec< kernel_t, impl_traits_t > Spec
Memory< p_assembled_t > assembled_data
void Assemble()
Partial assembly of quadrature point data.
void MultUnassembled(const Vector &x, Vector &y) const
void AssembleMatrix(SparseMatrix &M) const
Assemble the operator in a SparseMatrix.
double * GetData(int k)
Definition: densemat.hpp:1201
ShapeEvaluator< solFE_type, IR, real_t > solShapeEval
static const int qpts
kernel_t::template CoefficientEval< IR, coeff_t, impl_traits_t >::Type coeff_eval_t
void TestElementwiseExtractAssemble(const Vector &x, Vector &y) const
static const int dim
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
meshType::FE_type meshFE_type
MFEM_ALWAYS_INLINE void ElementAddMultAssembled(int el, solFieldEval &solFEval) const
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
static const int sdim
void MultAssembled(const Vector &x, Vector &y) const
static const int dofs
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
static const int SS
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
solFESpace::FE_type solFE_type
void AssembleFromSerializedNodes(const Vector &sNodes)
Partial assembly from "serialized" nodes.
static const int AB
void Serialize(const Vector &x, Vector &sx) const
complex_t - dof/qpt data type, real_t - ShapeEvaluator (FE basis) data type
Definition: tevaluator.hpp:997
TElementTransformation< meshType, IR, real_t > Trans_t
impl_traits_t::vcomplex_t vcomplex_t
double a
Definition: lissajous.cpp:41
void MakeDataOwner() const
Set the Vector data (host pointer) ownership flag.
Definition: vector.hpp:185
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
integ_t::coefficient_type coeff_t
coeff_t might be TConstantCoefficient or TFunctionCoefficient, for example
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...
TBilinearForm(const IntegratorType &integ, const FiniteElementSpace &sol_fes)
FieldEvaluator< solFESpace, solVecLayout_t, IR, complex_t, real_t > solFieldEval
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
int dim
Definition: ex24.cpp:53
void MultAssembledSerialized(const Vector &sx, Vector &sy) const
serialized vector sx –> serialized vector &#39;sy&#39;
void AddMult(DenseTensor &M, const Vector &x, Vector &y) const
Multiplication using assembled element matrices stored as a DenseTensor.
IntegratorType
Definition: par_example.cpp:39
Vector data type.
Definition: vector.hpp:58
static const int TE
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.
Definition: densemat.cpp:1591
RefCoord s[3]
solVecLayout_t solVecLayout
double * Data()
Definition: densemat.hpp:1213
impl_traits_t::vreal_t vreal_t
Spec::ElementMatrix ElementMatrix
data_t data[aligned_size >0?aligned_size:1]
Definition: ttensor.hpp:294
static const int BE
Abstract operator.
Definition: operator.hpp:24
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:1096
IntegratorType integ_t
impl_traits_t impl_traits_type
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...