MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tbilininteg.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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_BILININTEG
13 #define MFEM_TEMPLATE_BILININTEG
14 
15 #include "../config/tconfig.hpp"
16 #include "tcoefficient.hpp"
17 #include "tbilinearform.hpp"
18 
19 namespace mfem
20 {
21 
22 // Templated local bilinear form integrator kernels, cf. bilininteg.?pp
23 
24 // The Integrator class combines a kernel and a coefficient
25 
26 template <typename coeff_t, template<int,int,typename> class kernel_t>
28 {
29 public:
31 
32  template <int SDim, int Dim, typename complex_t>
33  struct kernel { typedef kernel_t<SDim,Dim,complex_t> type; };
34 
36 
37  TIntegrator(const coefficient_type &c) : coeff(c) { }
38 };
39 
40 
41 // Mass kernel
42 
43 template <int SDim, int Dim, typename complex_t>
45 {
46  typedef complex_t complex_type;
47 
48  // needed for the TElementTransformation::Result class
49  static const bool uses_Jacobians = true;
50 
51  // needed for the FieldEvaluator::Data class
52  static const bool in_values = true;
53  static const bool in_gradients = false;
54  static const bool out_values = true;
55  static const bool out_gradients = false;
56 
57  // Partially assembled data type for one element with the given number of
58  // quadrature points. This type is used in partial assembly, and partially
59  // assembled action.
60  template <int qpts>
62 
63  // Partially assembled data type for one element with the given number of
64  // quadrature points. This type is used in full element matrix assembly.
65  template <int qpts>
67 
68  template <typename IR, typename coeff_t, int NE>
70  {
72  };
73 
74  // Method used for un-assembled (matrix free) action.
75  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
76  // Q - CoefficientEval<>::Type
77  // q - CoefficientEval<>::Type::result_t
78  // val_qpts [M x NC x NE] - in/out data member in R
79  //
80  // val_qpts *= w det(J)
81  template <typename T_result_t, typename Q_t, typename q_t,
82  typename S_data_t>
83  static inline MFEM_ALWAYS_INLINE
84  void Action(const int k, const T_result_t &F,
85  const Q_t &Q, const q_t &q, S_data_t &R)
86  {
87  typedef typename T_result_t::Jt_type::data_type real_t;
88  const int M = S_data_t::eval_type::qpts;
89  const int NC = S_data_t::eval_type::vdim;
90  MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
91  "incompatible dimensions");
92  MFEM_FLOPS_ADD(M*(1+NC)); // TDet counts its flops
93  for (int i = 0; i < M; i++)
94  {
95  const complex_t wi =
96  Q.get(q,i,k) * TDet<real_t>(F.Jt.layout.ind14(i,k), F.Jt);
97  for (int j = 0; j < NC; j++)
98  {
99  R.val_qpts(i,j,k) *= wi;
100  }
101  }
102  }
103 
104  // Method defining partial assembly.
105  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
106  // Q - CoefficientEval<>::Type
107  // q - CoefficientEval<>::Type::result_t
108  // A [M] - partially assembled scalars
109  //
110  // A = w det(J)
111  template <typename T_result_t, typename Q_t, typename q_t, int qpts>
112  static inline MFEM_ALWAYS_INLINE
113  void Assemble(const int k, const T_result_t &F,
114  const Q_t &Q, const q_t &q, TVector<qpts,complex_t> &A)
115  {
116  typedef typename T_result_t::Jt_type::data_type real_t;
117 
118  const int M = T_result_t::Jt_type::layout_type::dim_1;
119  MFEM_STATIC_ASSERT(qpts == M, "incompatible dimensions");
120  MFEM_FLOPS_ADD(M); // TDet counts its flops
121  for (int i = 0; i < M; i++)
122  {
123  A[i] = Q.get(q,i,k) * TDet<real_t>(F.Jt.layout.ind14(i,k), F.Jt);
124  }
125  }
126 
127  // Method for partially assembled action.
128  // A [M] - partially assembled scalars
129  // val_qpts [M x NC x NE] - in/out data member in R
130  //
131  // val_qpts *= A
132  template <int qpts, typename S_data_t>
133  static inline MFEM_ALWAYS_INLINE
134  void MultAssembled(const int k, const TVector<qpts,complex_t> &A, S_data_t &R)
135  {
136  const int M = S_data_t::eval_type::qpts;
137  const int NC = S_data_t::eval_type::vdim;
138  MFEM_STATIC_ASSERT(qpts == M, "incompatible dimensions");
139  MFEM_FLOPS_ADD(M*NC);
140  for (int i = 0; i < M; i++)
141  {
142  for (int j = 0; j < NC; j++)
143  {
144  R.val_qpts(i,j,k) *= A[i];
145  }
146  }
147  }
148 };
149 
150 
151 // Diffusion kernel
152 
153 // complex_t - type for the assembled data
154 template <int SDim, int Dim, typename complex_t>
156 
157 // Diffusion kernel in 1D
158 template <typename complex_t>
159 struct TDiffusionKernel<1,1,complex_t>
160 {
161  typedef complex_t complex_type;
162 
163  // needed for the TElementTransformation::Result class
164  static const bool uses_Jacobians = true;
165 
166  // needed for the FieldEvaluator::Data class
167  static const bool in_values = false;
168  static const bool in_gradients = true;
169  static const bool out_values = false;
170  static const bool out_gradients = true;
171 
172  // Partially assembled data type for one element with the given number of
173  // quadrature points. This type is used in partial assembly, and partially
174  // assembled action.
175  template <int qpts>
176  struct p_asm_data { typedef TMatrix<qpts,1,complex_t> type; };
177 
178  // Partially assembled data type for one element with the given number of
179  // quadrature points. This type is used in full element matrix assembly.
180  template <int qpts>
181  struct f_asm_data { typedef TTensor3<qpts,1,1,complex_t> type; };
182 
183  template <typename IR, typename coeff_t, int NE>
184  struct CoefficientEval
185  {
187  };
188 
189  // Method used for un-assembled (matrix free) action.
190  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
191  // Q - CoefficientEval<>::Type
192  // q - CoefficientEval<>::Type::result_t
193  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
194  //
195  // grad_qpts = (w/det(J)) adj(J) adj(J)^t grad_qpts
196  template <typename T_result_t, typename Q_t, typename q_t,
197  typename S_data_t>
198  static inline MFEM_ALWAYS_INLINE
199  void Action(const int k, const T_result_t &F,
200  const Q_t &Q, const q_t &q, S_data_t &R)
201  {
202  const int M = S_data_t::eval_type::qpts;
203  const int NC = S_data_t::eval_type::vdim;
204  MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
205  "incompatible dimensions");
206  MFEM_FLOPS_ADD(M*(1+NC));
207  for (int i = 0; i < M; i++)
208  {
209  const complex_t wi = Q.get(q,i,k) / F.Jt(i,0,0,k);
210  for (int j = 0; j < NC; j++)
211  {
212  R.grad_qpts(i,0,j,k) *= wi;
213  }
214  }
215  }
216 
217  // Method defining partial assembly. The pointwise Dim x Dim matrices are
218  // stored as symmetric (when asm_type == p_asm_data, i.e. A.layout.rank == 2)
219  // or non-symmetric (when asm_type == f_asm_data, i.e. A.layout.rank == 3)
220  // matrices.
221  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
222  // Q - CoefficientEval<>::Type
223  // q - CoefficientEval<>::Type::result_t
224  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symm. matrices
225  // A [M x Dim x Dim] - partially assembled Dim x Dim matrices
226  //
227  // A = (w/det(J)) adj(J) adj(J)^t
228  template <typename T_result_t, typename Q_t, typename q_t, typename asm_type>
229  static inline MFEM_ALWAYS_INLINE
230  void Assemble(const int k, const T_result_t &F,
231  const Q_t &Q, const q_t &q, asm_type &A)
232  {
233  const int M = T_result_t::Jt_type::layout_type::dim_1;
234  MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
235  "incompatible dimensions");
236  MFEM_FLOPS_ADD(M);
237  for (int i = 0; i < M; i++)
238  {
239  // A[i] is A(i,0) or A(i,0,0)
240  A[i] = Q.get(q,i,k) / F.Jt(i,0,0,k);
241  }
242  }
243 
244  // Method for partially assembled action.
245  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symmetric
246  // matrices
247  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
248  //
249  // grad_qpts = A grad_qpts
250  template <int qpts, typename S_data_t>
251  static inline MFEM_ALWAYS_INLINE
252  void MultAssembled(const int k, const TMatrix<qpts,1,complex_t> &A,
253  S_data_t &R)
254  {
255  const int M = S_data_t::eval_type::qpts;
256  const int NC = S_data_t::eval_type::vdim;
257  MFEM_STATIC_ASSERT(qpts == M, "incompatible dimensions");
258  MFEM_FLOPS_ADD(M*NC);
259  for (int i = 0; i < M; i++)
260  {
261  for (int j = 0; j < NC; j++)
262  {
263  R.grad_qpts(i,0,j,k) *= A(i,0);
264  }
265  }
266  }
267 };
268 
269 // Diffusion kernel in 2D
270 template <typename complex_t>
271 struct TDiffusionKernel<2,2,complex_t>
272 {
273  typedef complex_t complex_type;
274 
275  // needed for the TElementTransformation::Result class
276  static const bool uses_Jacobians = true;
277 
278  // needed for the FieldEvaluator::Data class
279  static const bool in_values = false;
280  static const bool in_gradients = true;
281  static const bool out_values = false;
282  static const bool out_gradients = true;
283 
284  // Partially assembled data type for one element with the given number of
285  // quadrature points. This type is used in partial assembly, and partially
286  // assembled action. Stores one symmetric 2 x 2 matrix per point.
287  template <int qpts>
288  struct p_asm_data { typedef TMatrix<qpts,3,complex_t> type; };
289 
290  // Partially assembled data type for one element with the given number of
291  // quadrature points. This type is used in full element matrix assembly.
292  // Stores one general (non-symmetric) 2 x 2 matrix per point.
293  template <int qpts>
294  struct f_asm_data { typedef TTensor3<qpts,2,2,complex_t> type; };
295 
296  template <typename IR, typename coeff_t, int NE>
297  struct CoefficientEval
298  {
300  };
301 
302  // Method used for un-assembled (matrix free) action.
303  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
304  // Q - CoefficientEval<>::Type
305  // q - CoefficientEval<>::Type::result_t
306  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
307  //
308  // grad_qpts = (w/det(J)) adj(J) adj(J)^t grad_qpts
309  template <typename T_result_t, typename Q_t, typename q_t,
310  typename S_data_t>
311  static inline MFEM_ALWAYS_INLINE
312  void Action(const int k, const T_result_t &F,
313  const Q_t &Q, const q_t &q, S_data_t &R)
314  {
315  const int M = S_data_t::eval_type::qpts;
316  const int NC = S_data_t::eval_type::vdim;
317  MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
318  "incompatible dimensions");
319  MFEM_FLOPS_ADD(M*(4+NC*14));
320  for (int i = 0; i < M; i++)
321  {
322  typedef typename T_result_t::Jt_type::data_type real_t;
323  const real_t J11 = F.Jt(i,0,0,k);
324  const real_t J12 = F.Jt(i,1,0,k);
325  const real_t J21 = F.Jt(i,0,1,k);
326  const real_t J22 = F.Jt(i,1,1,k);
327  const complex_t w_det_J = Q.get(q,i,k) / (J11 * J22 - J21 * J12);
328  for (int j = 0; j < NC; j++)
329  {
330  const complex_t x1 = R.grad_qpts(i,0,j,k);
331  const complex_t x2 = R.grad_qpts(i,1,j,k);
332  // z = adj(J)^t x
333  const complex_t z1 = J22 * x1 - J21 * x2;
334  const complex_t z2 = J11 * x2 - J12 * x1;
335  R.grad_qpts(i,0,j,k) = w_det_J * (J22 * z1 - J12 * z2);
336  R.grad_qpts(i,1,j,k) = w_det_J * (J11 * z2 - J21 * z1);
337  }
338  }
339  }
340 
341  // Method defining partial assembly. The pointwise Dim x Dim matrices are
342  // stored as symmetric (when asm_type == p_asm_data, i.e. A.layout.rank == 2)
343  // or non-symmetric (when asm_type == f_asm_data, i.e. A.layout.rank == 3)
344  // matrices.
345  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
346  // Q - CoefficientEval<>::Type
347  // q - CoefficientEval<>::Type::result_t
348  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symm. matrices
349  // A [M x Dim x Dim] - partially assembled Dim x Dim matrices
350  //
351  // A = (w/det(J)) adj(J) adj(J)^t
352  template <typename T_result_t, typename Q_t, typename q_t, typename asm_type>
353  static inline MFEM_ALWAYS_INLINE
354  void Assemble(const int k, const T_result_t &F,
355  const Q_t &Q, const q_t &q, asm_type &A)
356  {
357  typedef typename T_result_t::Jt_type::data_type real_t;
358  const int M = T_result_t::Jt_type::layout_type::dim_1;
359  MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
360  "incompatible dimensions");
361  MFEM_FLOPS_ADD(16*M);
362  const bool Symm = (asm_type::layout_type::rank == 2);
363  for (int i = 0; i < M; i++)
364  {
365  const real_t J11 = F.Jt(i,0,0,k);
366  const real_t J12 = F.Jt(i,1,0,k);
367  const real_t J21 = F.Jt(i,0,1,k);
368  const real_t J22 = F.Jt(i,1,1,k);
369  const complex_t w_det_J = Q.get(q,i,k) / (J11 * J22 - J21 * J12);
370  internal::MatrixOps<2,2>::Symm<Symm>::Set(
371  A.layout.ind1(i), A,
372  + w_det_J * (J12*J12 + J22*J22), // (1,1)
373  - w_det_J * (J11*J12 + J21*J22), // (2,1)
374  + w_det_J * (J11*J11 + J21*J21) // (2,2)
375  );
376  }
377  }
378 
379  // Method for partially assembled action.
380  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symmetric
381  // matrices
382  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
383  //
384  // grad_qpts = A grad_qpts
385  template <int qpts, typename S_data_t>
386  static inline MFEM_ALWAYS_INLINE
387  void MultAssembled(const int k, const TMatrix<qpts,3,complex_t> &A,
388  S_data_t &R)
389  {
390  const int M = S_data_t::eval_type::qpts;
391  const int NC = S_data_t::eval_type::vdim;
392  MFEM_STATIC_ASSERT(qpts == M, "incompatible dimensions");
393  MFEM_FLOPS_ADD(6*M*NC);
394  for (int i = 0; i < M; i++)
395  {
396  const complex_t A11 = A(i,0);
397  const complex_t A21 = A(i,1);
398  const complex_t A22 = A(i,2);
399  for (int j = 0; j < NC; j++)
400  {
401  const complex_t x1 = R.grad_qpts(i,0,j,k);
402  const complex_t x2 = R.grad_qpts(i,1,j,k);
403  R.grad_qpts(i,0,j,k) = A11 * x1 + A21 * x2;
404  R.grad_qpts(i,1,j,k) = A21 * x1 + A22 * x2;
405  }
406  }
407  }
408 };
409 
410 // Diffusion kernel in 3D
411 template <typename complex_t>
412 struct TDiffusionKernel<3,3,complex_t>
413 {
414  typedef complex_t complex_type;
415 
416  // needed for the TElementTransformation::Result class
417  static const bool uses_Jacobians = true;
418 
419  // needed for the FieldEvaluator::Data class
420  static const bool in_values = false;
421  static const bool in_gradients = true;
422  static const bool out_values = false;
423  static const bool out_gradients = true;
424 
425  // Partially assembled data type for one element with the given number of
426  // quadrature points. This type is used in partial assembly, and partially
427  // assembled action. Stores one symmetric 3 x 3 matrix per point.
428  template <int qpts>
429  struct p_asm_data { typedef TMatrix<qpts,6,complex_t> type; };
430 
431  // Partially assembled data type for one element with the given number of
432  // quadrature points. This type is used in full element matrix assembly.
433  // Stores one general (non-symmetric) 3 x 3 matrix per point.
434  template <int qpts>
435  struct f_asm_data { typedef TTensor3<qpts,3,3,complex_t> type; };
436 
437  template <typename IR, typename coeff_t, int NE>
438  struct CoefficientEval
439  {
441  };
442 
443  // Method used for un-assembled (matrix free) action.
444  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
445  // Q - CoefficientEval<>::Type
446  // q - CoefficientEval<>::Type::result_t
447  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
448  //
449  // grad_qpts = (w/det(J)) adj(J) adj(J)^t grad_qpts
450  template <typename T_result_t, typename Q_t, typename q_t,
451  typename S_data_t>
452  static inline MFEM_ALWAYS_INLINE
453  void Action(const int k, const T_result_t &F,
454  const Q_t &Q, const q_t &q, S_data_t &R)
455  {
456  const int M = S_data_t::eval_type::qpts;
457  const int NC = S_data_t::eval_type::vdim;
458  MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
459  "incompatible dimensions");
460  MFEM_FLOPS_ADD(M); // just need to count Q/detJ
461  for (int i = 0; i < M; i++)
462  {
463  typedef typename T_result_t::Jt_type::data_type real_t;
464  TMatrix<3,3,real_t> adj_J;
465  const complex_t w_det_J =
466  (Q.get(q,i,k) /
467  TAdjDet<real_t>(F.Jt.layout.ind14(i,k).transpose_12(), F.Jt,
468  adj_J.layout, adj_J));
469  TMatrix<3,NC,complex_t> z; // z = adj(J)^t x
470  sMult_AB<false>(adj_J.layout.transpose_12(), adj_J,
471  R.grad_qpts.layout.ind14(i,k), R.grad_qpts,
472  z.layout, z);
473  z.Scale(w_det_J);
474  sMult_AB<false>(adj_J.layout, adj_J,
475  z.layout, z,
476  R.grad_qpts.layout.ind14(i,k), R.grad_qpts);
477  }
478  }
479 
480  // Method defining partial assembly. The pointwise Dim x Dim matrices are
481  // stored as symmetric (when asm_type == p_asm_data, i.e. A.layout.rank == 2)
482  // or non-symmetric (when asm_type == f_asm_data, i.e. A.layout.rank == 3)
483  // matrices.
484  // Jt [M x Dim x SDim x NE] - Jacobian transposed, data member in F
485  // Q - CoefficientEval<>::Type
486  // q - CoefficientEval<>::Type::result_t
487  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symm. matrices
488  // A [M x Dim x Dim] - partially assembled Dim x Dim matrices
489  //
490  // A = (w/det(J)) adj(J) adj(J)^t
491  template <typename T_result_t, typename Q_t, typename q_t, typename asm_type>
492  static inline MFEM_ALWAYS_INLINE
493  void Assemble(const int k, const T_result_t &F,
494  const Q_t &Q, const q_t &q, asm_type &A)
495  {
496  typedef typename T_result_t::Jt_type::data_type real_t;
497  const int M = T_result_t::Jt_type::layout_type::dim_1;
498  MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
499  "incompatible dimensions");
500  MFEM_FLOPS_ADD(37*M);
501  const bool Symm = (asm_type::layout_type::rank == 2);
502  for (int i = 0; i < M; i++)
503  {
504  TMatrix<3,3,real_t> B; // = adj(J)
505  const complex_t u =
506  (Q.get(q,i,k) /
507  TAdjDet<real_t>(F.Jt.layout.ind14(i,k).transpose_12(), F.Jt,
508  B.layout, B));
509  internal::MatrixOps<3,3>::Symm<Symm>::Set(
510  A.layout.ind1(i), A,
511  u*(B(0,0)*B(0,0)+B(0,1)*B(0,1)+B(0,2)*B(0,2)), // 1,1
512  u*(B(0,0)*B(1,0)+B(0,1)*B(1,1)+B(0,2)*B(1,2)), // 2,1
513  u*(B(0,0)*B(2,0)+B(0,1)*B(2,1)+B(0,2)*B(2,2)), // 3,1
514  u*(B(1,0)*B(1,0)+B(1,1)*B(1,1)+B(1,2)*B(1,2)), // 2,2
515  u*(B(1,0)*B(2,0)+B(1,1)*B(2,1)+B(1,2)*B(2,2)), // 3,2
516  u*(B(2,0)*B(2,0)+B(2,1)*B(2,1)+B(2,2)*B(2,2)) // 3,3
517  );
518  }
519  }
520 
521  // Method for partially assembled action.
522  // A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symmetric
523  // matrices
524  // grad_qpts [M x SDim x NC x NE] - in/out data member in R
525  //
526  // grad_qpts = A grad_qpts
527  template <int qpts, typename S_data_t>
528  static inline MFEM_ALWAYS_INLINE
529  void MultAssembled(const int k, const TMatrix<qpts,6,complex_t> &A,
530  S_data_t &R)
531  {
532  const int M = S_data_t::eval_type::qpts;
533  const int NC = S_data_t::eval_type::vdim;
534  MFEM_STATIC_ASSERT(qpts == M, "incompatible dimensions");
535  MFEM_FLOPS_ADD(15*M*NC);
536  for (int i = 0; i < M; i++)
537  {
538  const complex_t A11 = A(i,0);
539  const complex_t A21 = A(i,1);
540  const complex_t A31 = A(i,2);
541  const complex_t A22 = A(i,3);
542  const complex_t A32 = A(i,4);
543  const complex_t A33 = A(i,5);
544  for (int j = 0; j < NC; j++)
545  {
546  const complex_t x1 = R.grad_qpts(i,0,j,k);
547  const complex_t x2 = R.grad_qpts(i,1,j,k);
548  const complex_t x3 = R.grad_qpts(i,2,j,k);
549  R.grad_qpts(i,0,j,k) = A11*x1 + A21*x2 + A31*x3;
550  R.grad_qpts(i,1,j,k) = A21*x1 + A22*x2 + A32*x3;
551  R.grad_qpts(i,2,j,k) = A31*x1 + A32*x2 + A33*x3;
552  }
553  }
554  }
555 };
556 
557 } // namespace mfem
558 
559 #endif // MFEM_TEMPLATE_BILININTEG
static const bool uses_Jacobians
Definition: tbilininteg.hpp:49
static const bool in_values
Definition: tbilininteg.hpp:52
static const layout_type layout
Definition: ttensor.hpp:316
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
IntRuleCoefficient< IR, coeff_t, NE >::Type Type
static const bool out_gradients
Definition: tbilininteg.hpp:55
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 3, complex_t > &A, S_data_t &R)
IntRuleCoefficient< IR, coeff_t, NE >::Type Type
Definition: tbilininteg.hpp:71
TConstantCoefficient coeff_t
Definition: ex1.cpp:57
kernel_t< SDim, Dim, complex_t > type
Definition: tbilininteg.hpp:33
static const bool in_gradients
Definition: tbilininteg.hpp:53
TVector< qpts, complex_t > type
Definition: tbilininteg.hpp:61
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
Definition: tbilininteg.hpp:84
complex_t complex_type
Definition: tbilininteg.hpp:46
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 1, complex_t > &A, S_data_t &R)
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, TVector< qpts, complex_t > &A)
IntRuleCoefficient< IR, coeff_t, NE >::Type Type
IntRuleCoefficient< IR, coeff_t, NE >::Type Type
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
coeff_t coefficient_type
Definition: tbilininteg.hpp:30
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 6, complex_t > &A, S_data_t &R)
TVector< qpts, complex_t > type
Definition: tbilininteg.hpp:66
static const bool out_values
Definition: tbilininteg.hpp:54
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TVector< qpts, complex_t > &A, S_data_t &R)
void Scale(const data_t scale)
Definition: ttensor.hpp:297
TIntegrator(const coefficient_type &c)
Definition: tbilininteg.hpp:37