MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tbilinearform.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_TEMPLATE_BILINEAR_FORM
13 #define MFEM_TEMPLATE_BILINEAR_FORM
14 
15 #include "../config/tconfig.hpp"
16 #include "../linalg/ttensor.hpp"
17 #include "bilinearform.hpp"
18 #include "tevaluator.hpp"
19 #include "teltrans.hpp"
20 #include "tcoefficient.hpp"
21 #include "fespace.hpp"
22 
23 namespace mfem
24 {
25 
26 // Templated bilinear form class, cf. bilinearform.?pp
27 
28 // complex_t - sol dof data type
29 // real_t - mesh nodes, sol basis, mesh basis data type
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>
34 class TBilinearForm : public Operator
35 {
36 protected:
37  typedef complex_t complex_type;
38  typedef real_t real_type;
39 
40  typedef typename meshType::FE_type meshFE_type;
42  typedef typename solFESpace::FE_type solFE_type;
44  typedef solVecLayout_t solVecLayout_type;
45 
46  static const int dim = meshType::dim;
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;
51 
52  typedef IntegratorType integ_t;
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;
57 
59  template <int NE> struct T_result
60  {
61  static const int EvalOps =
62  Trans_t::template Get<coeff_t,kernel_t>::EvalOps;
63  typedef typename Trans_t::template Result<EvalOps,NE> Type;
64  };
65 
66  typedef FieldEvaluator<solFESpace,solVecLayout_t,IR,
67  complex_t,real_t> solFieldEval;
68  template <int BE> struct S_spec
69  {
70  typedef typename solFieldEval::template Spec<kernel_t,BE> Spec;
71  typedef typename Spec::DataType DataType;
72  typedef typename Spec::ElementMatrix ElementMatrix;
73  };
74 
75  // Data members
76 
77  meshType mesh;
79 
82  mutable solFESpace solFES;
83  solVecLayout_t solVecLayout;
84 
86 
88 
90 
91 public:
92  TBilinearForm(const IntegratorType &integ, const FiniteElementSpace &sol_fes)
93  : Operator(sol_fes.GetNDofs()*vdim),
94  mesh(*sol_fes.GetMesh()),
95  meshEval(mesh.fe),
96  sol_fe(*sol_fes.FEColl()),
97  solEval(sol_fe),
98  solFES(sol_fe, sol_fes),
99  solVecLayout(sol_fes),
100  int_rule(),
101  coeff(integ.coeff),
102  assembled_data(NULL)
103  { }
104 
105  virtual ~TBilinearForm()
106  {
107  delete [] assembled_data;
108  }
109 
110  virtual void Mult(const Vector &x, Vector &y) const
111  {
112  if (assembled_data)
113  {
114  const int num_elem = 1;
115  MultAssembled<num_elem>(x, y);
116  }
117  else
118  {
119  MultUnassembled(x, y);
120  }
121  }
122 
123  // complex_t = double
124  void MultUnassembled(const Vector &x, Vector &y) const
125  {
126  y = 0.0;
127 
128  const int BE = 1; // batch-size of elements
129  typedef typename kernel_t::template
130  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
131 
132  // For better performance, create stack copies of solFES, and solEval
133  // inside 'solFEval'. The element-transformation 'T' also copies the
134  // meshFES, meshEval, etc internally.
135  // Is performance actually better with this implementation?
136  Trans_t T(mesh, meshEval);
138  x.GetData(), y.GetData());
139  coeff_eval_t wQ(int_rule, coeff);
140 
141  const int NE = mesh.GetNE();
142  for (int el = 0; el < NE; el++)
143  {
144 #if 0
145  typename S_spec<BE>::DataType R;
146  solFEval.Eval(el, R);
147 
148  typename T_result<BE>::Type F;
149  T.Eval(el, F);
150 #else
151  typename T_result<BE>::Type F;
152  T.Eval(el, F);
153 
154  typename S_spec<BE>::DataType R;
155  solFEval.Eval(el, R);
156 #endif
157 
158  typename coeff_eval_t::result_t res;
159  wQ.Eval(F, res);
160 
161  kernel_t::Action(0, F, wQ, res, R);
162 
163  solFEval.template Assemble<true>(R);
164  }
165  }
166 
167  // Partial assembly of quadrature point data
168  void Assemble()
169  {
170  const int BE = 1; // batch-size of elements
171  typedef typename kernel_t::template
172  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
173 
174  Trans_t T(mesh, meshEval);
175  coeff_eval_t wQ(int_rule, coeff);
176 
177  const int NE = mesh.GetNE();
178  if (!assembled_data)
179  {
180  assembled_data = new p_assembled_t[NE];
181  }
182  for (int el = 0; el < NE; el++) // BE == 1
183  {
184  typename T_result<BE>::Type F;
185  T.Eval(el, F);
186 
187  typename coeff_eval_t::result_t res;
188  wQ.Eval(F, res);
189 
190  for (int k = 0; k < BE; k++)
191  {
192  kernel_t::Assemble(k, F, wQ, res, assembled_data[el+k]);
193  }
194  }
195  }
196 
197  template <int num_elem>
198  inline MFEM_ALWAYS_INLINE
199  void ElementAddMultAssembled(int el, solFieldEval &solFEval) const
200  {
201  typename S_spec<num_elem>::DataType R;
202  solFEval.Eval(el, R);
203 
204  for (int k = 0; k < num_elem; k++)
205  {
206  kernel_t::MultAssembled(k, assembled_data[el+k], R);
207  }
208 
209  solFEval.template Assemble<true>(R);
210  }
211 
212  // complex_t = double
213  template <int num_elem>
214  void MultAssembled(const Vector &x, Vector &y) const
215  {
216  y = 0.0;
217 
219  x.GetData(), y.GetData());
220 
221  const int NE = mesh.GetNE();
222  const int bNE = NE-NE%num_elem;
223  for (int el = 0; el < bNE; el += num_elem)
224  {
225  ElementAddMultAssembled<num_elem>(el, solFEval);
226  }
227  for (int el = bNE; el < NE; el++)
228  {
229  ElementAddMultAssembled<1>(el, solFEval);
230  }
231  }
232 
233 #ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
234  // complex_t = double
236  {
237  y = 0.0;
238 
240  solFESpace solFES(this->solFES);
241 
243 
244  const int NE = mesh.GetNE();
245  for (int el = 0; el < NE; el++)
246  {
247  solFES.SetElement(el);
248 
249  solFES.VectorExtract(solVecLayout, x, xy_dof.layout, xy_dof);
250  solFES.VectorAssemble(xy_dof.layout, xy_dof, solVecLayout, y);
251  }
252  }
253 
254  // real_t = double
255  void SerializeNodes(Vector &sNodes) const
256  {
257  typedef typename meshType::FESpace_type meshFESpace;
258  meshFESpace meshFES(mesh.t_fes);
260 
261  const int NE = mesh.GetNE();
262  sNodes.SetSize(lnodes_t::size*NE);
263  real_t *lNodes = sNodes.GetData();
264  for (int el = 0; el < NE; el++)
265  {
266  meshFES.SetElement(el);
267  meshFES.VectorExtract(mesh.node_layout, mesh.Nodes,
268  lnodes_t::layout, lNodes);
269  lNodes += lnodes_t::size;
270  }
271  }
272 
273  // partial assembly from "serialized" nodes
274  // real_t = double
276  {
277  const int BE = 1; // batch-size of elements
278  typedef typename kernel_t::template
279  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
280 
281  Trans_t T(this->mesh, this->meshEval);
282  coeff_eval_t wQ(int_rule, coeff);
283 
284  const int NE = mesh.GetNE();
285  if (!assembled_data)
286  {
287  assembled_data = new p_assembled_t[NE];
288  }
289  for (int el = 0; el < NE; el++)
290  {
291  typename T_result<BE>::Type F;
292  T.EvalSerialized(el, sNodes.GetData(), F);
293 
294  typename coeff_eval_t::result_t res;
295  wQ.Eval(F, res);
296 
297  kernel_t::Assemble(0, F, wQ, res, assembled_data[el]);
298  }
299  }
300 
301  // complex_t = double
302  void Serialize(const Vector &x, Vector &sx) const
303  {
304  solVecLayout_t solVecLayout(this->solVecLayout);
305  typedef TTensor3<dofs,vdim,1,complex_t> vdof_data_t;
306  solFESpace solFES(this->solFES);
307 
308  const int NE = mesh.GetNE();
309  sx.SetSize(vdim*dofs*NE);
310  complex_t *loc_sx = sx.GetData();
311  for (int el = 0; el < NE; el++)
312  {
313  solFES.SetElement(el);
314  solFES.VectorExtract(solVecLayout, x, vdof_data_t::layout, loc_sx);
315  loc_sx += vdim*dofs;
316  }
317  }
318 
319  // serialized vector sx --> serialized vector 'sy'
320  // complex_t = double
321  void MultAssembledSerialized(const Vector &sx, Vector &sy) const
322  {
323  solFieldEval solFEval(solFES, solEval, solVecLayout, NULL, NULL);
324 
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++)
329  {
330  typename S_spec<1>::DataType R;
331  solFEval.EvalSerialized(loc_sx, R);
332 
333  kernel_t::MultAssembled(0, assembled_data[el], R);
334 
335  solFEval.template AssembleSerialized<false>(R, loc_sy);
336 
337  loc_sx += vdim*dofs;
338  loc_sy += vdim*dofs;
339  }
340  }
341 #endif // MFEM_TEMPLATE_ENABLE_SERIALIZE
342 
343  // Assemble the operator in a SparseMatrix.
344  // complex_t = double
346  {
347  const int BE = 1; // batch-size of elements
348  typedef typename kernel_t::template
349  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
350 
351  Trans_t T(mesh, meshEval);
352  solFESpace solFES(this->solFES);
354  solVecLayout_t solVecLayout(this->solVecLayout);
355  coeff_eval_t wQ(int_rule, coeff);
356 
357  const int NE = mesh.GetNE();
358  for (int el = 0; el < NE; el++)
359  {
360  f_assembled_t asm_qpt_data;
361  {
362  typename T_result<BE>::Type F;
363  T.Eval(el, F);
364 
365  typename coeff_eval_t::result_t res;
366  wQ.Eval(F, res);
367 
368  kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
369  }
370 
371  // For now, when vdim > 1, assume block-diagonal matrix with the same
372  // diagonal block for all components.
373  TMatrix<dofs,dofs> M_loc;
375  asm_qpt_data.layout, asm_qpt_data, M_loc.layout, M_loc, solEval);
376 
377  solFES.SetElement(el);
378  for (int bi = 0; bi < vdim; bi++)
379  {
380  solFES.AssembleBlock(bi, bi, solVecLayout, M_loc, M);
381  }
382  }
383  }
384 
385  // Assemble element matrices and store them as a DenseTensor object.
386  // complex_t = double
388  {
389  const int BE = 1; // batch-size of elements
390  typedef typename kernel_t::template
391  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
392 
393  Trans_t T(mesh, meshEval);
395  coeff_eval_t wQ(int_rule, coeff);
396 
397  const int NE = mesh.GetNE();
398  for (int el = 0; el < NE; el++)
399  {
400  f_assembled_t asm_qpt_data;
401  {
402  typename T_result<BE>::Type F;
403  T.Eval(el, F);
404 
405  typename coeff_eval_t::result_t res;
406  wQ.Eval(F, res);
407 
408  kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
409  }
410 
411  // For now, when vdim > 1, assume block-diagonal matrix with the same
412  // diagonal block for all components.
413  // M is assumed to be (dof x dof x NE).
414  TMatrix<dofs,dofs> M_loc;
416  asm_qpt_data.layout, asm_qpt_data, M_loc.layout, M_loc, solEval);
417 
418  complex_t *M_data = M.GetData(el);
419  M_loc.template AssignTo<AssignOp::Set>(M_data);
420  }
421  }
422 
423  // Assemble element matrices and add them to the bilinear form
424  // complex_t = double
426  {
427  const int BE = 1; // batch-size of elements
428  typedef typename kernel_t::template
429  CoefficientEval<IR,coeff_t,BE>::Type coeff_eval_t;
430 
431  Trans_t T(mesh, meshEval);
433  coeff_eval_t wQ(int_rule, coeff);
434 
435  Array<int> vdofs;
436  const Array<int> *dof_map = sol_fe.GetDofMap();
437  const int *dof_map_ = dof_map->GetData();
438  DenseMatrix M_loc_perm(dofs*vdim,dofs*vdim); // initialized with zeros
439 
440  const int NE = mesh.GetNE();
441  for (int el = 0; el < NE; el++)
442  {
443  f_assembled_t asm_qpt_data;
444  {
445  typename T_result<BE>::Type F;
446  T.Eval(el, F);
447 
448  typename coeff_eval_t::result_t res;
449  wQ.Eval(F, res);
450 
451  kernel_t::Assemble(0, F, wQ, res, asm_qpt_data);
452  }
453 
454  // For now, when vdim > 1, assume block-diagonal matrix with the same
455  // diagonal block for all components.
456  TMatrix<dofs,dofs> M_loc;
458  asm_qpt_data.layout, asm_qpt_data, M_loc.layout, M_loc, solEval);
459 
460  if (dof_map) // switch from tensor-product ordering
461  {
462  for (int i = 0; i < dofs; i++)
463  {
464  for (int j = 0; j < dofs; j++)
465  {
466  M_loc_perm(dof_map_[i],dof_map_[j]) = M_loc(i,j);
467  }
468  }
469  for (int bi = 1; bi < vdim; bi++)
470  {
471  M_loc_perm.CopyMN(M_loc_perm, dofs, dofs, 0, 0,
472  bi*dofs, bi*dofs);
473  }
474  a.AssembleElementMatrix(el, M_loc_perm, vdofs);
475  }
476  else
477  {
478  DenseMatrix DM(M_loc.data, dofs, dofs);
479  if (vdim == 1)
480  {
481  a.AssembleElementMatrix(el, DM, vdofs);
482  }
483  else
484  {
485  for (int bi = 0; bi < vdim; bi++)
486  {
487  M_loc_perm.CopyMN(DM, dofs, dofs, 0, 0, bi*dofs, bi*dofs);
488  }
489  a.AssembleElementMatrix(el, M_loc_perm, vdofs);
490  }
491  }
492  }
493  }
494 
495  // Multiplication using assembled element matrices stored as a DenseTensor.
496  // complex_t = double
497  void AddMult(DenseTensor &M, const Vector &x, Vector &y) const
498  {
499  // For now, when vdim > 1, assume block-diagonal matrix with the same
500  // diagonal block for all components.
501  // M is assumed to be (dof x dof x NE).
502  solVecLayout_t solVecLayout(this->solVecLayout);
503  const int NE = mesh.GetNE();
504  for (int el = 0; el < NE; el++)
505  {
506  TTensor3<dofs,vdim,1,complex_t> x_dof, y_dof;
507 
508  solFES.SetElement(el);
509  solFES.VectorExtract(solVecLayout, x, x_dof.layout, x_dof);
510  Mult_AB<false>(TMatrix<dofs,dofs>::layout,
511  M(el).Data(),
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);
515  }
516  }
517 };
518 
519 } // namespace mfem
520 
521 #endif // MFEM_TEMPLATE_BILINEAR_FORM
meshType::FE_type meshFE_type
FieldEvaluator< solFESpace, solVecLayout_t, IR, complex_t, real_t > solFieldEval
float real_t
Definition: mesh.cpp:4276
MFEM_ALWAYS_INLINE void EvalSerialized(int el, const real_t *nodeData, Result< EvalOps, NE > &F)
Definition: teltrans.hpp:119
void SetSize(int s)
Resize the vector if the new size is different.
Definition: vector.hpp:263
static const layout_type layout
Definition: ttensor.hpp:311
solVecLayout_t solVecLayout
MFEM_ALWAYS_INLINE void Eval(DataType &F)
T * GetData()
Returns the data.
Definition: array.hpp:91
integ_t::coefficient_type coeff_t
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
void AssembleBilinearForm(BilinearForm &a) const
kernel_t::template p_asm_data< qpts >::type p_assembled_t
void AddMult(DenseTensor &M, const Vector &x, Vector &y) const
static const layout_type layout
Definition: ttensor.hpp:347
double * GetData() const
Definition: vector.hpp:90
void AssembleFromSerializedNodes(const Vector &sNodes)
static const int dim
void MultUnassembled(const Vector &x, Vector &y) const
int dim
Definition: ex3.cpp:47
TBilinearForm(const IntegratorType &integ, const FiniteElementSpace &sol_fes)
Data type sparse matrix.
Definition: sparsemat.hpp:38
MFEM_ALWAYS_INLINE void EvalSerialized(const complex_t *loc_dofs, DataType &F)
solFESpace::FE_type solFE_type
double * GetData(int k)
Definition: densemat.hpp:624
void Serialize(const Vector &x, Vector &sx) const
Trans_t::template Result< EvalOps, NE > Type
meshShapeEval meshEval
solVecLayout_t solVecLayout_type
ShapeEvaluator< solFE_type, IR, real_t > solShapeEval
void AssembleElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
void TestElementwiseExtractAssemble(const Vector &x, Vector &y) const
static const int dofs
MFEM_ALWAYS_INLINE void ElementAddMultAssembled(int el, solFieldEval &solFEval) const
static const int sdim
MFEM_ALWAYS_INLINE void Eval(int el, Result< EvalOps, NE > &F)
Definition: teltrans.hpp:111
IntegratorType integ_t
TElementTransformation< meshType, IR, real_t > Trans_t
void MultAssembledSerialized(const Vector &sx, Vector &sy) const
solFieldEval::template Spec< kernel_t, BE > Spec
void AssembleMatrix(DenseTensor &M) const
integ_t::template kernel< sdim, dim, complex_t >::type kernel_t
static const int qpts
Vector data type.
Definition: vector.hpp:33
ShapeEvaluator< meshFE_type, IR, real_t > meshShapeEval
void SerializeNodes(Vector &sNodes) const
static const int vdim
p_assembled_t * assembled_data
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:2460
double * Data()
Definition: densemat.hpp:626
void MultAssembled(const Vector &x, Vector &y) const
data_t data[aligned_size >0?aligned_size:1]
Definition: ttensor.hpp:249
solShapeEval solEval
kernel_t::template f_asm_data< qpts >::type f_assembled_t
Abstract operator.
Definition: operator.hpp:21
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:569
void AssembleMatrix(SparseMatrix &M) const
Spec::ElementMatrix ElementMatrix