MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_rt.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 // Raviart-Thomas Finite Element classes
13 
14 #include "fe_rt.hpp"
15 #include "../coefficient.hpp"
16 
17 namespace mfem
18 {
19 
20 using namespace std;
21 
22 const double RT_QuadrilateralElement::nk[8] =
23 { 0., -1., 1., 0., 0., 1., -1., 0. };
24 
26  const int cb_type,
27  const int ob_type)
28  : VectorTensorFiniteElement(2, 2*(p + 1)*(p + 2), p + 1, cb_type, ob_type,
29  H_DIV, DofMapType::L2_DOF_MAP),
30  dof2nk(dof),
31  cp(poly1d.ClosedPoints(p + 1, cb_type))
32 {
33  if (obasis1d.IsIntegratedType()) { is_nodal = false; }
34 
36 
37  const double *op = poly1d.OpenPoints(p, ob_type);
38  const int dof2 = dof/2;
39 
40 #ifndef MFEM_THREAD_SAFE
41  shape_cx.SetSize(p + 2);
42  shape_ox.SetSize(p + 1);
43  shape_cy.SetSize(p + 2);
44  shape_oy.SetSize(p + 1);
45  dshape_cx.SetSize(p + 2);
46  dshape_cy.SetSize(p + 2);
47 #endif
48 
49  // edges
50  int o = 0;
51  for (int i = 0; i <= p; i++) // (0,1)
52  {
53  dof_map[1*dof2 + i + 0*(p + 1)] = o++;
54  }
55  for (int i = 0; i <= p; i++) // (1,2)
56  {
57  dof_map[0*dof2 + (p + 1) + i*(p + 2)] = o++;
58  }
59  for (int i = 0; i <= p; i++) // (2,3)
60  {
61  dof_map[1*dof2 + (p - i) + (p + 1)*(p + 1)] = o++;
62  }
63  for (int i = 0; i <= p; i++) // (3,0)
64  {
65  dof_map[0*dof2 + 0 + (p - i)*(p + 2)] = o++;
66  }
67 
68  // interior
69  for (int j = 0; j <= p; j++) // x-components
70  for (int i = 1; i <= p; i++)
71  {
72  dof_map[0*dof2 + i + j*(p + 2)] = o++;
73  }
74  for (int j = 1; j <= p; j++) // y-components
75  for (int i = 0; i <= p; i++)
76  {
77  dof_map[1*dof2 + i + j*(p + 1)] = o++;
78  }
79 
80  // dof orientations
81  // x-components
82  for (int j = 0; j <= p; j++)
83  for (int i = 0; i <= p/2; i++)
84  {
85  int idx = 0*dof2 + i + j*(p + 2);
86  dof_map[idx] = -1 - dof_map[idx];
87  }
88  if (p%2 == 1)
89  for (int j = p/2 + 1; j <= p; j++)
90  {
91  int idx = 0*dof2 + (p/2 + 1) + j*(p + 2);
92  dof_map[idx] = -1 - dof_map[idx];
93  }
94  // y-components
95  for (int j = 0; j <= p/2; j++)
96  for (int i = 0; i <= p; i++)
97  {
98  int idx = 1*dof2 + i + j*(p + 1);
99  dof_map[idx] = -1 - dof_map[idx];
100  }
101  if (p%2 == 1)
102  for (int i = 0; i <= p/2; i++)
103  {
104  int idx = 1*dof2 + i + (p/2 + 1)*(p + 1);
105  dof_map[idx] = -1 - dof_map[idx];
106  }
107 
108  o = 0;
109  for (int j = 0; j <= p; j++)
110  for (int i = 0; i <= p + 1; i++)
111  {
112  int idx;
113  if ((idx = dof_map[o++]) < 0)
114  {
115  idx = -1 - idx;
116  dof2nk[idx] = 3;
117  }
118  else
119  {
120  dof2nk[idx] = 1;
121  }
122  Nodes.IntPoint(idx).Set2(cp[i], op[j]);
123  }
124  for (int j = 0; j <= p + 1; j++)
125  for (int i = 0; i <= p; i++)
126  {
127  int idx;
128  if ((idx = dof_map[o++]) < 0)
129  {
130  idx = -1 - idx;
131  dof2nk[idx] = 0;
132  }
133  else
134  {
135  dof2nk[idx] = 2;
136  }
137  Nodes.IntPoint(idx).Set2(op[i], cp[j]);
138  }
139 }
140 
142  DenseMatrix &shape) const
143 {
144  const int pp1 = order;
145 
146 #ifdef MFEM_THREAD_SAFE
147  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
148  Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
149 #endif
150 
152  {
153  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
154  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
155  obasis1d.ScaleIntegrated(false);
156  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
157  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
158  }
159  else
160  {
161  cbasis1d.Eval(ip.x, shape_cx);
162  cbasis1d.Eval(ip.y, shape_cy);
163  obasis1d.Eval(ip.x, shape_ox);
164  obasis1d.Eval(ip.y, shape_oy);
165  }
166 
167  int o = 0;
168  for (int j = 0; j < pp1; j++)
169  for (int i = 0; i <= pp1; i++)
170  {
171  int idx, s;
172  if ((idx = dof_map[o++]) < 0)
173  {
174  idx = -1 - idx, s = -1;
175  }
176  else
177  {
178  s = +1;
179  }
180  shape(idx,0) = s*shape_cx(i)*shape_oy(j);
181  shape(idx,1) = 0.;
182  }
183  for (int j = 0; j <= pp1; j++)
184  for (int i = 0; i < pp1; i++)
185  {
186  int idx, s;
187  if ((idx = dof_map[o++]) < 0)
188  {
189  idx = -1 - idx, s = -1;
190  }
191  else
192  {
193  s = +1;
194  }
195  shape(idx,0) = 0.;
196  shape(idx,1) = s*shape_ox(i)*shape_cy(j);
197  }
198 }
199 
201  Vector &divshape) const
202 {
203  const int pp1 = order;
204 
205 #ifdef MFEM_THREAD_SAFE
206  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
207  Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
208 #endif
209 
210  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
211  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
213  {
214  obasis1d.ScaleIntegrated(false);
215  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
216  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
217  }
218  else
219  {
220  obasis1d.Eval(ip.x, shape_ox);
221  obasis1d.Eval(ip.y, shape_oy);
222  }
223 
224  int o = 0;
225  for (int j = 0; j < pp1; j++)
226  for (int i = 0; i <= pp1; i++)
227  {
228  int idx, s;
229  if ((idx = dof_map[o++]) < 0)
230  {
231  idx = -1 - idx, s = -1;
232  }
233  else
234  {
235  s = +1;
236  }
237  divshape(idx) = s*dshape_cx(i)*shape_oy(j);
238  }
239  for (int j = 0; j <= pp1; j++)
240  for (int i = 0; i < pp1; i++)
241  {
242  int idx, s;
243  if ((idx = dof_map[o++]) < 0)
244  {
245  idx = -1 - idx, s = -1;
246  }
247  else
248  {
249  s = +1;
250  }
251  divshape(idx) = s*shape_ox(i)*dshape_cy(j);
252  }
253 }
254 
257  Vector &dofs) const
258 {
259  MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
260  double vk[Geometry::MaxDim];
261  Vector xk(vk, vc.GetVDim());
262 
264  const int nqpt = ir.GetNPoints();
265 
266  IntegrationPoint ip2d;
267 
268  int o = 0;
269  for (int c = 0; c < 2; c++)
270  {
271  int im = (c == 0) ? order + 1 : order;
272  int jm = (c == 1) ? order + 1 : order;
273  for (int j = 0; j < jm; j++)
274  for (int i = 0; i < im; i++)
275  {
276  int idx = dof_map[o++];
277  if (idx < 0) { idx = -1 - idx; }
278  int ic = (c == 0) ? j : i;
279  const double h = cp[ic+1] - cp[ic];
280  double val = 0.0;
281  for (int k = 0; k < nqpt; k++)
282  {
283  const IntegrationPoint &ip1d = ir.IntPoint(k);
284  if (c == 0) { ip2d.Set2(cp[i], cp[j] + (h*ip1d.x)); }
285  else { ip2d.Set2(cp[i] + (h*ip1d.x), cp[j]); }
286  Trans.SetIntPoint(&ip2d);
287  vc.Eval(xk, Trans, ip2d);
288  // nk^t adj(J) xk
289  const double ipval = Trans.AdjugateJacobian().InnerProduct(vk,
290  nk + dof2nk[idx]*dim);
291  val += ip1d.weight*ipval;
292  }
293  dofs(idx) = val*h;
294  }
295  }
296 }
297 
298 
299 const double RT_HexahedronElement::nk[18] =
300 { 0.,0.,-1., 0.,-1.,0., 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,0.,1. };
301 
303  const int cb_type,
304  const int ob_type)
305  : VectorTensorFiniteElement(3, 3*(p + 1)*(p + 1)*(p + 2), p + 1, cb_type,
306  ob_type, H_DIV, DofMapType::L2_DOF_MAP),
307  dof2nk(dof),
308  cp(poly1d.ClosedPoints(p + 1, cb_type))
309 {
310  if (obasis1d.IsIntegratedType()) { is_nodal = false; }
311 
313 
314  const double *op = poly1d.OpenPoints(p, ob_type);
315  const int dof3 = dof/3;
316 
317 #ifndef MFEM_THREAD_SAFE
318  shape_cx.SetSize(p + 2);
319  shape_ox.SetSize(p + 1);
320  shape_cy.SetSize(p + 2);
321  shape_oy.SetSize(p + 1);
322  shape_cz.SetSize(p + 2);
323  shape_oz.SetSize(p + 1);
324  dshape_cx.SetSize(p + 2);
325  dshape_cy.SetSize(p + 2);
326  dshape_cz.SetSize(p + 2);
327 #endif
328 
329  // faces
330  int o = 0;
331  for (int j = 0; j <= p; j++) // (3,2,1,0) -- bottom
332  for (int i = 0; i <= p; i++)
333  {
334  dof_map[2*dof3 + i + ((p - j) + 0*(p + 1))*(p + 1)] = o++;
335  }
336  for (int j = 0; j <= p; j++) // (0,1,5,4) -- front
337  for (int i = 0; i <= p; i++)
338  {
339  dof_map[1*dof3 + i + (0 + j*(p + 2))*(p + 1)] = o++;
340  }
341  for (int j = 0; j <= p; j++) // (1,2,6,5) -- right
342  for (int i = 0; i <= p; i++)
343  {
344  dof_map[0*dof3 + (p + 1) + (i + j*(p + 1))*(p + 2)] = o++;
345  }
346  for (int j = 0; j <= p; j++) // (2,3,7,6) -- back
347  for (int i = 0; i <= p; i++)
348  {
349  dof_map[1*dof3 + (p - i) + ((p + 1) + j*(p + 2))*(p + 1)] = o++;
350  }
351  for (int j = 0; j <= p; j++) // (3,0,4,7) -- left
352  for (int i = 0; i <= p; i++)
353  {
354  dof_map[0*dof3 + 0 + ((p - i) + j*(p + 1))*(p + 2)] = o++;
355  }
356  for (int j = 0; j <= p; j++) // (4,5,6,7) -- top
357  for (int i = 0; i <= p; i++)
358  {
359  dof_map[2*dof3 + i + (j + (p + 1)*(p + 1))*(p + 1)] = o++;
360  }
361 
362  // interior
363  // x-components
364  for (int k = 0; k <= p; k++)
365  for (int j = 0; j <= p; j++)
366  for (int i = 1; i <= p; i++)
367  {
368  dof_map[0*dof3 + i + (j + k*(p + 1))*(p + 2)] = o++;
369  }
370  // y-components
371  for (int k = 0; k <= p; k++)
372  for (int j = 1; j <= p; j++)
373  for (int i = 0; i <= p; i++)
374  {
375  dof_map[1*dof3 + i + (j + k*(p + 2))*(p + 1)] = o++;
376  }
377  // z-components
378  for (int k = 1; k <= p; k++)
379  for (int j = 0; j <= p; j++)
380  for (int i = 0; i <= p; i++)
381  {
382  dof_map[2*dof3 + i + (j + k*(p + 1))*(p + 1)] = o++;
383  }
384 
385  // dof orientations
386  // for odd p, do not change the orientations in the mid-planes
387  // {i = p/2 + 1}, {j = p/2 + 1}, {k = p/2 + 1} in the x, y, z-components
388  // respectively.
389  // x-components
390  for (int k = 0; k <= p; k++)
391  for (int j = 0; j <= p; j++)
392  for (int i = 0; i <= p/2; i++)
393  {
394  int idx = 0*dof3 + i + (j + k*(p + 1))*(p + 2);
395  dof_map[idx] = -1 - dof_map[idx];
396  }
397  // y-components
398  for (int k = 0; k <= p; k++)
399  for (int j = 0; j <= p/2; j++)
400  for (int i = 0; i <= p; i++)
401  {
402  int idx = 1*dof3 + i + (j + k*(p + 2))*(p + 1);
403  dof_map[idx] = -1 - dof_map[idx];
404  }
405  // z-components
406  for (int k = 0; k <= p/2; k++)
407  for (int j = 0; j <= p; j++)
408  for (int i = 0; i <= p; i++)
409  {
410  int idx = 2*dof3 + i + (j + k*(p + 1))*(p + 1);
411  dof_map[idx] = -1 - dof_map[idx];
412  }
413 
414  o = 0;
415  // x-components
416  for (int k = 0; k <= p; k++)
417  for (int j = 0; j <= p; j++)
418  for (int i = 0; i <= p + 1; i++)
419  {
420  int idx;
421  if ((idx = dof_map[o++]) < 0)
422  {
423  idx = -1 - idx;
424  dof2nk[idx] = 4;
425  }
426  else
427  {
428  dof2nk[idx] = 2;
429  }
430  Nodes.IntPoint(idx).Set3(cp[i], op[j], op[k]);
431  }
432  // y-components
433  for (int k = 0; k <= p; k++)
434  for (int j = 0; j <= p + 1; j++)
435  for (int i = 0; i <= p; i++)
436  {
437  int idx;
438  if ((idx = dof_map[o++]) < 0)
439  {
440  idx = -1 - idx;
441  dof2nk[idx] = 1;
442  }
443  else
444  {
445  dof2nk[idx] = 3;
446  }
447  Nodes.IntPoint(idx).Set3(op[i], cp[j], op[k]);
448  }
449  // z-components
450  for (int k = 0; k <= p + 1; k++)
451  for (int j = 0; j <= p; j++)
452  for (int i = 0; i <= p; i++)
453  {
454  int idx;
455  if ((idx = dof_map[o++]) < 0)
456  {
457  idx = -1 - idx;
458  dof2nk[idx] = 0;
459  }
460  else
461  {
462  dof2nk[idx] = 5;
463  }
464  Nodes.IntPoint(idx).Set3(op[i], op[j], cp[k]);
465  }
466 }
467 
469  DenseMatrix &shape) const
470 {
471  const int pp1 = order;
472 
473 #ifdef MFEM_THREAD_SAFE
474  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
475  Vector shape_cz(pp1 + 1), shape_oz(pp1);
476  Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
477 #endif
478 
480  {
481  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
482  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
483  cbasis1d.Eval(ip.z, shape_cz, dshape_cz);
484  obasis1d.ScaleIntegrated(false);
485  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
486  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
487  obasis1d.EvalIntegrated(dshape_cz, shape_oz);
488  }
489  else
490  {
491  cbasis1d.Eval(ip.x, shape_cx);
492  cbasis1d.Eval(ip.y, shape_cy);
493  cbasis1d.Eval(ip.z, shape_cz);
494  obasis1d.Eval(ip.x, shape_ox);
495  obasis1d.Eval(ip.y, shape_oy);
496  obasis1d.Eval(ip.z, shape_oz);
497  }
498 
499  int o = 0;
500  // x-components
501  for (int k = 0; k < pp1; k++)
502  for (int j = 0; j < pp1; j++)
503  for (int i = 0; i <= pp1; i++)
504  {
505  int idx, s;
506  if ((idx = dof_map[o++]) < 0)
507  {
508  idx = -1 - idx, s = -1;
509  }
510  else
511  {
512  s = +1;
513  }
514  shape(idx,0) = s*shape_cx(i)*shape_oy(j)*shape_oz(k);
515  shape(idx,1) = 0.;
516  shape(idx,2) = 0.;
517  }
518  // y-components
519  for (int k = 0; k < pp1; k++)
520  for (int j = 0; j <= pp1; j++)
521  for (int i = 0; i < pp1; i++)
522  {
523  int idx, s;
524  if ((idx = dof_map[o++]) < 0)
525  {
526  idx = -1 - idx, s = -1;
527  }
528  else
529  {
530  s = +1;
531  }
532  shape(idx,0) = 0.;
533  shape(idx,1) = s*shape_ox(i)*shape_cy(j)*shape_oz(k);
534  shape(idx,2) = 0.;
535  }
536  // z-components
537  for (int k = 0; k <= pp1; k++)
538  for (int j = 0; j < pp1; j++)
539  for (int i = 0; i < pp1; i++)
540  {
541  int idx, s;
542  if ((idx = dof_map[o++]) < 0)
543  {
544  idx = -1 - idx, s = -1;
545  }
546  else
547  {
548  s = +1;
549  }
550  shape(idx,0) = 0.;
551  shape(idx,1) = 0.;
552  shape(idx,2) = s*shape_ox(i)*shape_oy(j)*shape_cz(k);
553  }
554 }
555 
557  Vector &divshape) const
558 {
559  const int pp1 = order;
560 
561 #ifdef MFEM_THREAD_SAFE
562  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
563  Vector shape_cz(pp1 + 1), shape_oz(pp1);
564  Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
565 #endif
566 
567  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
568  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
569  cbasis1d.Eval(ip.z, shape_cz, dshape_cz);
571  {
572  obasis1d.ScaleIntegrated(false);
573  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
574  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
575  obasis1d.EvalIntegrated(dshape_cz, shape_oz);
576  }
577  else
578  {
579  obasis1d.Eval(ip.x, shape_ox);
580  obasis1d.Eval(ip.y, shape_oy);
581  obasis1d.Eval(ip.z, shape_oz);
582  }
583 
584  int o = 0;
585  // x-components
586  for (int k = 0; k < pp1; k++)
587  for (int j = 0; j < pp1; j++)
588  for (int i = 0; i <= pp1; i++)
589  {
590  int idx, s;
591  if ((idx = dof_map[o++]) < 0)
592  {
593  idx = -1 - idx, s = -1;
594  }
595  else
596  {
597  s = +1;
598  }
599  divshape(idx) = s*dshape_cx(i)*shape_oy(j)*shape_oz(k);
600  }
601  // y-components
602  for (int k = 0; k < pp1; k++)
603  for (int j = 0; j <= pp1; j++)
604  for (int i = 0; i < pp1; i++)
605  {
606  int idx, s;
607  if ((idx = dof_map[o++]) < 0)
608  {
609  idx = -1 - idx, s = -1;
610  }
611  else
612  {
613  s = +1;
614  }
615  divshape(idx) = s*shape_ox(i)*dshape_cy(j)*shape_oz(k);
616  }
617  // z-components
618  for (int k = 0; k <= pp1; k++)
619  for (int j = 0; j < pp1; j++)
620  for (int i = 0; i < pp1; i++)
621  {
622  int idx, s;
623  if ((idx = dof_map[o++]) < 0)
624  {
625  idx = -1 - idx, s = -1;
626  }
627  else
628  {
629  s = +1;
630  }
631  divshape(idx) = s*shape_ox(i)*shape_oy(j)*dshape_cz(k);
632  }
633 }
634 
637  Vector &dofs) const
638 {
639  MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
640  double vq[Geometry::MaxDim];
641  Vector xq(vq, vc.GetVDim());
642 
644  const int nqpt = ir2d.GetNPoints();
645 
646  IntegrationPoint ip3d;
647 
648  int o = 0;
649  for (int c = 0; c < 3; c++)
650  {
651  int im = (c == 0) ? order + 1 : order;
652  int jm = (c == 1) ? order + 1 : order;
653  int km = (c == 2) ? order + 1 : order;
654  for (int k = 0; k < km; k++)
655  for (int j = 0; j < jm; j++)
656  for (int i = 0; i < im; i++)
657  {
658  int idx = dof_map[o++];
659  if (idx < 0) { idx = -1 - idx; }
660  int ic1, ic2;
661  if (c == 0) { ic1 = j; ic2 = k; }
662  else if (c == 1) { ic1 = i; ic2 = k; }
663  else { ic1 = i; ic2 = j; }
664  const double h1 = cp[ic1+1] - cp[ic1];
665  const double h2 = cp[ic2+1] - cp[ic2];
666  double val = 0.0;
667  for (int q = 0; q < nqpt; q++)
668  {
669  const IntegrationPoint &ip2d = ir2d.IntPoint(q);
670  if (c == 0) { ip3d.Set3(cp[i], cp[j] + h1*ip2d.x, cp[k] + h2*ip2d.y); }
671  else if (c == 1) { ip3d.Set3(cp[i] + h1*ip2d.x, cp[j], cp[k] + h2*ip2d.y); }
672  else { ip3d.Set3(cp[i] + h1*ip2d.x, cp[j] + h2*ip2d.y, cp[k]); }
673  Trans.SetIntPoint(&ip3d);
674  vc.Eval(xq, Trans, ip3d);
675  // nk^t adj(J) xq
676  const double ipval
677  = Trans.AdjugateJacobian().InnerProduct(vq, nk + dof2nk[idx]*dim);
678  val += ip2d.weight*ipval;
679  }
680  dofs(idx) = val*h1*h2;
681  }
682  }
683 }
684 
685 
686 const double RT_TriangleElement::nk[6] =
687 { 0., -1., 1., 1., -1., 0. };
688 
689 const double RT_TriangleElement::c = 1./3.;
690 
692  : VectorFiniteElement(2, Geometry::TRIANGLE, (p + 1)*(p + 3), p + 1,
693  H_DIV, FunctionSpace::Pk),
694  dof2nk(dof)
695 {
696  const double *iop = (p > 0) ? poly1d.OpenPoints(p - 1) : NULL;
697  const double *bop = poly1d.OpenPoints(p);
698 
699 #ifndef MFEM_THREAD_SAFE
700  shape_x.SetSize(p + 1);
701  shape_y.SetSize(p + 1);
702  shape_l.SetSize(p + 1);
703  dshape_x.SetSize(p + 1);
704  dshape_y.SetSize(p + 1);
705  dshape_l.SetSize(p + 1);
706  u.SetSize(dof, dim);
707  divu.SetSize(dof);
708 #else
709  Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
710 #endif
711 
712  // edges
713  int o = 0;
714  for (int i = 0; i <= p; i++) // (0,1)
715  {
716  Nodes.IntPoint(o).Set2(bop[i], 0.);
717  dof2nk[o++] = 0;
718  }
719  for (int i = 0; i <= p; i++) // (1,2)
720  {
721  Nodes.IntPoint(o).Set2(bop[p-i], bop[i]);
722  dof2nk[o++] = 1;
723  }
724  for (int i = 0; i <= p; i++) // (2,0)
725  {
726  Nodes.IntPoint(o).Set2(0., bop[p-i]);
727  dof2nk[o++] = 2;
728  }
729 
730  // interior
731  for (int j = 0; j < p; j++)
732  for (int i = 0; i + j < p; i++)
733  {
734  double w = iop[i] + iop[j] + iop[p-1-i-j];
735  Nodes.IntPoint(o).Set2(iop[i]/w, iop[j]/w);
736  dof2nk[o++] = 0;
737  Nodes.IntPoint(o).Set2(iop[i]/w, iop[j]/w);
738  dof2nk[o++] = 2;
739  }
740 
741  DenseMatrix T(dof);
742  for (int k = 0; k < dof; k++)
743  {
744  const IntegrationPoint &ip = Nodes.IntPoint(k);
745  poly1d.CalcBasis(p, ip.x, shape_x);
746  poly1d.CalcBasis(p, ip.y, shape_y);
747  poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l);
748  const double *n_k = nk + 2*dof2nk[k];
749 
750  o = 0;
751  for (int j = 0; j <= p; j++)
752  for (int i = 0; i + j <= p; i++)
753  {
754  double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
755  T(o++, k) = s*n_k[0];
756  T(o++, k) = s*n_k[1];
757  }
758  for (int i = 0; i <= p; i++)
759  {
760  double s = shape_x(i)*shape_y(p-i);
761  T(o++, k) = s*((ip.x - c)*n_k[0] + (ip.y - c)*n_k[1]);
762  }
763  }
764 
765  Ti.Factor(T);
766  // mfem::out << "RT_TriangleElement(" << p << ") : "; Ti.TestInversion();
767 }
768 
770  DenseMatrix &shape) const
771 {
772  const int p = order - 1;
773 
774 #ifdef MFEM_THREAD_SAFE
775  Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
776  DenseMatrix u(dof, dim);
777 #endif
778 
779  poly1d.CalcBasis(p, ip.x, shape_x);
780  poly1d.CalcBasis(p, ip.y, shape_y);
781  poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l);
782 
783  int o = 0;
784  for (int j = 0; j <= p; j++)
785  for (int i = 0; i + j <= p; i++)
786  {
787  double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
788  u(o,0) = s; u(o,1) = 0; o++;
789  u(o,0) = 0; u(o,1) = s; o++;
790  }
791  for (int i = 0; i <= p; i++)
792  {
793  double s = shape_x(i)*shape_y(p-i);
794  u(o,0) = (ip.x - c)*s;
795  u(o,1) = (ip.y - c)*s;
796  o++;
797  }
798 
799  Ti.Mult(u, shape);
800 }
801 
803  Vector &divshape) const
804 {
805  const int p = order - 1;
806 
807 #ifdef MFEM_THREAD_SAFE
808  Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
809  Vector dshape_x(p + 1), dshape_y(p + 1), dshape_l(p + 1);
810  Vector divu(dof);
811 #endif
812 
813  poly1d.CalcBasis(p, ip.x, shape_x, dshape_x);
814  poly1d.CalcBasis(p, ip.y, shape_y, dshape_y);
815  poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l, dshape_l);
816 
817  int o = 0;
818  for (int j = 0; j <= p; j++)
819  for (int i = 0; i + j <= p; i++)
820  {
821  int k = p - i - j;
822  divu(o++) = (dshape_x(i)*shape_l(k) -
823  shape_x(i)*dshape_l(k))*shape_y(j);
824  divu(o++) = (dshape_y(j)*shape_l(k) -
825  shape_y(j)*dshape_l(k))*shape_x(i);
826  }
827  for (int i = 0; i <= p; i++)
828  {
829  int j = p - i;
830  divu(o++) = ((shape_x(i) + (ip.x - c)*dshape_x(i))*shape_y(j) +
831  (shape_y(j) + (ip.y - c)*dshape_y(j))*shape_x(i));
832  }
833 
834  Ti.Mult(divu, divshape);
835 }
836 
837 
838 const double RT_TetrahedronElement::nk[12] =
839 { 1,1,1, -1,0,0, 0,-1,0, 0,0,-1 };
840 // { .5,.5,.5, -.5,0,0, 0,-.5,0, 0,0,-.5}; // n_F |F|
841 
842 const double RT_TetrahedronElement::c = 1./4.;
843 
845  : VectorFiniteElement(3, Geometry::TETRAHEDRON, (p + 1)*(p + 2)*(p + 4)/2,
846  p + 1, H_DIV, FunctionSpace::Pk),
847  dof2nk(dof)
848 {
849  const double *iop = (p > 0) ? poly1d.OpenPoints(p - 1) : NULL;
850  const double *bop = poly1d.OpenPoints(p);
851 
852 #ifndef MFEM_THREAD_SAFE
853  shape_x.SetSize(p + 1);
854  shape_y.SetSize(p + 1);
855  shape_z.SetSize(p + 1);
856  shape_l.SetSize(p + 1);
857  dshape_x.SetSize(p + 1);
858  dshape_y.SetSize(p + 1);
859  dshape_z.SetSize(p + 1);
860  dshape_l.SetSize(p + 1);
861  u.SetSize(dof, dim);
862  divu.SetSize(dof);
863 #else
864  Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
865 #endif
866 
867  int o = 0;
868  // faces (see Mesh::GenerateFaces in mesh/mesh.cpp,
869  // the constructor of H1_TetrahedronElement)
870  for (int j = 0; j <= p; j++)
871  for (int i = 0; i + j <= p; i++) // (1,2,3)
872  {
873  double w = bop[i] + bop[j] + bop[p-i-j];
874  Nodes.IntPoint(o).Set3(bop[p-i-j]/w, bop[i]/w, bop[j]/w);
875  dof2nk[o++] = 0;
876  }
877  for (int j = 0; j <= p; j++)
878  for (int i = 0; i + j <= p; i++) // (0,3,2)
879  {
880  double w = bop[i] + bop[j] + bop[p-i-j];
881  Nodes.IntPoint(o).Set3(0., bop[j]/w, bop[i]/w);
882  dof2nk[o++] = 1;
883  }
884  for (int j = 0; j <= p; j++)
885  for (int i = 0; i + j <= p; i++) // (0,1,3)
886  {
887  double w = bop[i] + bop[j] + bop[p-i-j];
888  Nodes.IntPoint(o).Set3(bop[i]/w, 0., bop[j]/w);
889  dof2nk[o++] = 2;
890  }
891  for (int j = 0; j <= p; j++)
892  for (int i = 0; i + j <= p; i++) // (0,2,1)
893  {
894  double w = bop[i] + bop[j] + bop[p-i-j];
895  Nodes.IntPoint(o).Set3(bop[j]/w, bop[i]/w, 0.);
896  dof2nk[o++] = 3;
897  }
898 
899  // interior
900  for (int k = 0; k < p; k++)
901  for (int j = 0; j + k < p; j++)
902  for (int i = 0; i + j + k < p; i++)
903  {
904  double w = iop[i] + iop[j] + iop[k] + iop[p-1-i-j-k];
905  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
906  dof2nk[o++] = 1;
907  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
908  dof2nk[o++] = 2;
909  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
910  dof2nk[o++] = 3;
911  }
912 
913  DenseMatrix T(dof);
914  for (int m = 0; m < dof; m++)
915  {
916  const IntegrationPoint &ip = Nodes.IntPoint(m);
917  poly1d.CalcBasis(p, ip.x, shape_x);
918  poly1d.CalcBasis(p, ip.y, shape_y);
919  poly1d.CalcBasis(p, ip.z, shape_z);
920  poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l);
921  const double *nm = nk + 3*dof2nk[m];
922 
923  o = 0;
924  for (int k = 0; k <= p; k++)
925  for (int j = 0; j + k <= p; j++)
926  for (int i = 0; i + j + k <= p; i++)
927  {
928  double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
929  T(o++, m) = s * nm[0];
930  T(o++, m) = s * nm[1];
931  T(o++, m) = s * nm[2];
932  }
933  for (int j = 0; j <= p; j++)
934  for (int i = 0; i + j <= p; i++)
935  {
936  double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
937  T(o++, m) = s*((ip.x - c)*nm[0] + (ip.y - c)*nm[1] +
938  (ip.z - c)*nm[2]);
939  }
940  }
941 
942  Ti.Factor(T);
943  // mfem::out << "RT_TetrahedronElement(" << p << ") : "; Ti.TestInversion();
944 }
945 
947  DenseMatrix &shape) const
948 {
949  const int p = order - 1;
950 
951 #ifdef MFEM_THREAD_SAFE
952  Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
953  DenseMatrix u(dof, dim);
954 #endif
955 
956  poly1d.CalcBasis(p, ip.x, shape_x);
957  poly1d.CalcBasis(p, ip.y, shape_y);
958  poly1d.CalcBasis(p, ip.z, shape_z);
959  poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l);
960 
961  int o = 0;
962  for (int k = 0; k <= p; k++)
963  for (int j = 0; j + k <= p; j++)
964  for (int i = 0; i + j + k <= p; i++)
965  {
966  double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
967  u(o,0) = s; u(o,1) = 0; u(o,2) = 0; o++;
968  u(o,0) = 0; u(o,1) = s; u(o,2) = 0; o++;
969  u(o,0) = 0; u(o,1) = 0; u(o,2) = s; o++;
970  }
971  for (int j = 0; j <= p; j++)
972  for (int i = 0; i + j <= p; i++)
973  {
974  double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
975  u(o,0) = (ip.x - c)*s; u(o,1) = (ip.y - c)*s; u(o,2) = (ip.z - c)*s;
976  o++;
977  }
978 
979  Ti.Mult(u, shape);
980 }
981 
983  Vector &divshape) const
984 {
985  const int p = order - 1;
986 
987 #ifdef MFEM_THREAD_SAFE
988  Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
989  Vector dshape_x(p + 1), dshape_y(p + 1), dshape_z(p + 1), dshape_l(p + 1);
990  Vector divu(dof);
991 #endif
992 
993  poly1d.CalcBasis(p, ip.x, shape_x, dshape_x);
994  poly1d.CalcBasis(p, ip.y, shape_y, dshape_y);
995  poly1d.CalcBasis(p, ip.z, shape_z, dshape_z);
996  poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l, dshape_l);
997 
998  int o = 0;
999  for (int k = 0; k <= p; k++)
1000  for (int j = 0; j + k <= p; j++)
1001  for (int i = 0; i + j + k <= p; i++)
1002  {
1003  int l = p - i - j - k;
1004  divu(o++) = (dshape_x(i)*shape_l(l) -
1005  shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
1006  divu(o++) = (dshape_y(j)*shape_l(l) -
1007  shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
1008  divu(o++) = (dshape_z(k)*shape_l(l) -
1009  shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
1010  }
1011  for (int j = 0; j <= p; j++)
1012  for (int i = 0; i + j <= p; i++)
1013  {
1014  int k = p - i - j;
1015  divu(o++) =
1016  (shape_x(i) + (ip.x - c)*dshape_x(i))*shape_y(j)*shape_z(k) +
1017  (shape_y(j) + (ip.y - c)*dshape_y(j))*shape_x(i)*shape_z(k) +
1018  (shape_z(k) + (ip.z - c)*dshape_z(k))*shape_x(i)*shape_y(j);
1019  }
1020 
1021  Ti.Mult(divu, divshape);
1022 }
1023 
1024 const double RT_WedgeElement::nk[15] =
1025 { 0,0,-1, 0,0,1, 0,-1,0, 1,1,0, -1,0,0};
1026 
1028  : VectorFiniteElement(3, Geometry::PRISM,
1029  (p + 2) * ((p + 1) * (p + 2)) / 2 +
1030  (p + 1) * (p + 1) * (p + 3), p + 1,
1031  H_DIV, FunctionSpace::Qk),
1032  dof2nk(dof),
1033  t_dof(dof),
1034  s_dof(dof),
1035  L2TriangleFE(p),
1036  RTTriangleFE(p),
1037  H1SegmentFE(p + 1),
1038  L2SegmentFE(p)
1039 {
1040  MFEM_ASSERT(L2TriangleFE.GetDof() * H1SegmentFE.GetDof() +
1041  RTTriangleFE.GetDof() * L2SegmentFE.GetDof() == dof,
1042  "Mismatch in number of degrees of freedom "
1043  "when building RT_WedgeElement!");
1044 
1045  const int pm1 = p - 1;
1046 
1047 #ifndef MFEM_THREAD_SAFE
1048  tl2_shape.SetSize(L2TriangleFE.GetDof());
1049  sh1_shape.SetSize(H1SegmentFE.GetDof());
1050  trt_shape.SetSize(RTTriangleFE.GetDof(), 2);
1051  sl2_shape.SetSize(L2SegmentFE.GetDof());
1052  sh1_dshape.SetSize(H1SegmentFE.GetDof(), 1);
1053  trt_dshape.SetSize(RTTriangleFE.GetDof());
1054 #endif
1055 
1056  const IntegrationRule &tl2_n = L2TriangleFE.GetNodes();
1057  const IntegrationRule &trt_n = RTTriangleFE.GetNodes();
1058  const IntegrationRule &sh1_n = H1SegmentFE.GetNodes();
1059  const IntegrationRule &sl2_n = L2SegmentFE.GetNodes();
1060 
1061  // faces
1062  int o = 0;
1063  int l = 0;
1064  // (0,2,1) -- bottom
1065  for (int j = 0; j <= p; j++)
1066  for (int i = 0; i + j <= p; i++)
1067  {
1068  l = j + i * (2 * p + 3 - i) / 2;
1069  t_dof[o] = l; s_dof[o] = 0; dof2nk[o] = 0;
1070  const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1071  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1072  o++;
1073  }
1074  // (3,4,5) -- top
1075  l = 0;
1076  for (int j = 0; j <= p; j++)
1077  for (int i = 0; i + j <= p; i++)
1078  {
1079  t_dof[o] = l; s_dof[o] = 1; dof2nk[o] = 1; l++;
1080  const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1081  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1082  o++;
1083  }
1084  // (0, 1, 4, 3) -- xz plane
1085  for (int j = 0; j <= p; j++)
1086  for (int i = 0; i <= p; i++)
1087  {
1088  t_dof[o] = i; s_dof[o] = j; dof2nk[o] = 2;
1089  const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1090  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1091  o++;
1092  }
1093  // (1, 2, 5, 4) -- (y-x)z plane
1094  for (int j = 0; j <= p; j++)
1095  for (int i = 0; i <= p; i++)
1096  {
1097  t_dof[o] = p + 1 + i; s_dof[o] = j; dof2nk[o] = 3;
1098  const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1099  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1100  o++;
1101  }
1102  // (2, 0, 3, 5) -- yz plane
1103  for (int j = 0; j <= p; j++)
1104  for (int i = 0; i <= p; i++)
1105  {
1106  t_dof[o] = 2 * p + 2 + i; s_dof[o] = j; dof2nk[o] = 4;
1107  const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1108  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1109  o++;
1110  }
1111 
1112  // interior
1113  for (int k = 0; k < L2SegmentFE.GetDof(); k++)
1114  {
1115  l = 0;
1116  for (int j = 0; j <= pm1; j++)
1117  for (int i = 0; i + j <= pm1; i++)
1118  {
1119  t_dof[o] = 3 * (p + 1) + 2 * l; s_dof[o] = k;
1120  dof2nk[o] = 2;
1121  const IntegrationPoint & t_ip0 = trt_n.IntPoint(t_dof[o]);
1122  const IntegrationPoint & s_ip0 = sl2_n.IntPoint(s_dof[o]);
1123  Nodes.IntPoint(o).Set3(t_ip0.x, t_ip0.y, s_ip0.x);
1124  o++;
1125  t_dof[o] = 3 * (p + 1) + 2 * l + 1; s_dof[o] = k;
1126  dof2nk[o] = 4; l++;
1127  const IntegrationPoint & t_ip1 = trt_n.IntPoint(t_dof[o]);
1128  const IntegrationPoint & s_ip1 = sl2_n.IntPoint(s_dof[o]);
1129  Nodes.IntPoint(o).Set3(t_ip1.x, t_ip1.y, s_ip1.x);
1130  o++;
1131  }
1132  }
1133  for (int k = 2; k < H1SegmentFE.GetDof(); k++)
1134  {
1135  for (l = 0; l < L2TriangleFE.GetDof(); l++)
1136  {
1137  t_dof[o] = l; s_dof[o] = k; dof2nk[o] = 1;
1138  const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1139  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1140  o++;
1141  }
1142  }
1143 }
1144 
1146  DenseMatrix &shape) const
1147 {
1148 #ifdef MFEM_THREAD_SAFE
1149  DenseMatrix trt_shape(RTTriangleFE.GetDof(), 2);
1150  Vector tl2_shape(L2TriangleFE.GetDof());
1151  Vector sh1_shape(H1SegmentFE.GetDof());
1152  Vector sl2_shape(L2SegmentFE.GetDof());
1153 #endif
1154 
1155  IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1156 
1157  L2TriangleFE.CalcShape(ip, tl2_shape);
1158  RTTriangleFE.CalcVShape(ip, trt_shape);
1159  H1SegmentFE.CalcShape(ipz, sh1_shape);
1160  L2SegmentFE.CalcShape(ipz, sl2_shape);
1161 
1162  for (int i=0; i<dof; i++)
1163  {
1164  if ( dof2nk[i] >= 2 )
1165  {
1166  shape(i, 0) = trt_shape(t_dof[i], 0) * sl2_shape[s_dof[i]];
1167  shape(i, 1) = trt_shape(t_dof[i], 1) * sl2_shape[s_dof[i]];
1168  shape(i, 2) = 0.0;
1169  }
1170  else
1171  {
1172  double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1173  shape(i, 0) = 0.0;
1174  shape(i, 1) = 0.0;
1175  shape(i, 2) = s * tl2_shape[t_dof[i]] * sh1_shape(s_dof[i]);
1176  }
1177  }
1178 }
1179 
1181  Vector &divshape) const
1182 {
1183 #ifdef MFEM_THREAD_SAFE
1184  Vector trt_dshape(RTTriangleFE.GetDof());
1185  Vector tl2_shape(L2TriangleFE.GetDof());
1186  Vector sl2_shape(L2SegmentFE.GetDof());
1187  DenseMatrix sh1_dshape(H1SegmentFE.GetDof(), 1);
1188 #endif
1189 
1190  IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1191 
1192  RTTriangleFE.CalcDivShape(ip, trt_dshape);
1193  L2TriangleFE.CalcShape(ip, tl2_shape);
1194 
1195  L2SegmentFE.CalcShape(ipz, sl2_shape);
1196  H1SegmentFE.CalcDShape(ipz, sh1_dshape);
1197 
1198  for (int i=0; i<dof; i++)
1199  {
1200  if ( dof2nk[i] >= 2 )
1201  {
1202  divshape(i) = trt_dshape(t_dof[i]) * sl2_shape(s_dof[i]);
1203  }
1204  else
1205  {
1206  double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1207  divshape(i) = s * tl2_shape(t_dof[i]) * sh1_dshape(s_dof[i], 0);
1208  }
1209  }
1210 }
1211 
1212 const double RT_R1D_SegmentElement::nk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1213 
1215  const int cb_type,
1216  const int ob_type)
1217  : VectorFiniteElement(1, Geometry::SEGMENT, 3 * p + 4, p + 1,
1218  H_DIV, FunctionSpace::Pk),
1219  dof2nk(dof),
1220  cbasis1d(poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1221  obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
1222 {
1223  // Override default dimension for VectorFiniteElements
1224  vdim = 3;
1225 
1226  const double *cp = poly1d.ClosedPoints(p + 1, cb_type);
1227  const double *op = poly1d.OpenPoints(p, ob_type);
1228 
1229 #ifndef MFEM_THREAD_SAFE
1230  shape_cx.SetSize(p + 2);
1231  shape_ox.SetSize(p + 1);
1232  dshape_cx.SetSize(p + 2);
1233 #endif
1234 
1235  dof_map.SetSize(dof);
1236 
1237  int o = 0;
1238  // nodes
1239  // (0)
1240  Nodes.IntPoint(o).x = cp[0]; // x-directed
1241  dof_map[0] = o; dof2nk[o++] = 0;
1242 
1243  // (1)
1244  Nodes.IntPoint(o).x = cp[p+1]; // x-directed
1245  dof_map[p+1] = o; dof2nk[o++] = 0;
1246 
1247  // interior
1248  // x-components
1249  for (int i = 1; i <= p; i++)
1250  {
1251  Nodes.IntPoint(o).x = cp[i];
1252  dof_map[i] = o; dof2nk[o++] = 0;
1253  }
1254  // y-components
1255  for (int i = 0; i <= p; i++)
1256  {
1257  Nodes.IntPoint(o).x = op[i];
1258  dof_map[p+i+2] = o; dof2nk[o++] = 1;
1259  }
1260  // z-components
1261  for (int i = 0; i <= p; i++)
1262  {
1263  Nodes.IntPoint(o).x = op[i];
1264  dof_map[2*p+3+i] = o; dof2nk[o++] = 2;
1265  }
1266 }
1267 
1269  DenseMatrix &shape) const
1270 {
1271  const int p = order;
1272 
1273 #ifdef MFEM_THREAD_SAFE
1274  Vector shape_cx(p + 1), shape_ox(p);
1275 #endif
1276 
1277  cbasis1d.Eval(ip.x, shape_cx);
1278  obasis1d.Eval(ip.x, shape_ox);
1279 
1280  int o = 0;
1281  // x-components
1282  for (int i = 0; i <= p; i++)
1283  {
1284  int idx = dof_map[o++];
1285  shape(idx,0) = shape_cx(i);
1286  shape(idx,1) = 0.;
1287  shape(idx,2) = 0.;
1288  }
1289  // y-components
1290  for (int i = 0; i < p; i++)
1291  {
1292  int idx = dof_map[o++];
1293  shape(idx,0) = 0.;
1294  shape(idx,1) = shape_ox(i);
1295  shape(idx,2) = 0.;
1296  }
1297  // z-components
1298  for (int i = 0; i < p; i++)
1299  {
1300  int idx = dof_map[o++];
1301  shape(idx,0) = 0.;
1302  shape(idx,1) = 0.;
1303  shape(idx,2) = shape_ox(i);
1304  }
1305 }
1306 
1308  DenseMatrix &shape) const
1309 {
1310  CalcVShape(Trans.GetIntPoint(), shape);
1311  const DenseMatrix & J = Trans.Jacobian();
1312  MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1313  "RT_R1D_SegmentElement cannot be embedded in "
1314  "2 or 3 dimensional spaces");
1315  for (int i=0; i<dof; i++)
1316  {
1317  shape(i, 0) *= J(0,0);
1318  }
1319  shape *= (1.0 / Trans.Weight());
1320 }
1321 
1323  Vector &divshape) const
1324 {
1325  const int p = order;
1326 
1327 #ifdef MFEM_THREAD_SAFE
1328  Vector shape_cx(p + 1);
1329  Vector dshape_cx(p + 1);
1330 #endif
1331 
1332  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
1333 
1334  int o = 0;
1335  // x-components
1336  for (int i = 0; i <= p; i++)
1337  {
1338  int idx = dof_map[o++];
1339  divshape(idx) = dshape_cx(i);
1340  }
1341  // y-components
1342  for (int i = 0; i < p; i++)
1343  {
1344  int idx = dof_map[o++];
1345  divshape(idx) = 0.;
1346  }
1347  // z-components
1348  for (int i = 0; i < p; i++)
1349  {
1350  int idx = dof_map[o++];
1351  divshape(idx) = 0.;
1352  }
1353 }
1354 
1357  Vector &dofs) const
1358 {
1359  double data[3];
1360  Vector vk1(data, 1);
1361  Vector vk3(data, 3);
1362 
1363  double * nk_ptr = const_cast<double*>(nk);
1364 
1365  for (int k = 0; k < dof; k++)
1366  {
1367  Trans.SetIntPoint(&Nodes.IntPoint(k));
1368 
1369  vc.Eval(vk3, Trans, Nodes.IntPoint(k));
1370  // dof_k = nk^t adj(J) vk
1371  Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1372  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1373 
1374  dofs(k) = Trans.AdjugateJacobian().InnerProduct(vk1, n1) +
1375  Trans.Weight() * vk3(1) * n3(1) +
1376  Trans.Weight() * vk3(2) * n3(2);
1377  }
1378 }
1379 
1382  DenseMatrix &I) const
1383 {
1384  if (fe.GetRangeType() == SCALAR)
1385  {
1386  double vk[Geometry::MaxDim];
1387  Vector shape(fe.GetDof());
1388 
1389  double * nk_ptr = const_cast<double*>(nk);
1390 
1391  I.SetSize(dof, vdim*fe.GetDof());
1392  for (int k = 0; k < dof; k++)
1393  {
1394  const IntegrationPoint &ip = Nodes.IntPoint(k);
1395 
1396  Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1397  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1398 
1399  fe.CalcShape(ip, shape);
1400  Trans.SetIntPoint(&ip);
1401  // Transform RT face normals from reference to physical space
1402  // vk = adj(J)^T nk
1403  Trans.AdjugateJacobian().MultTranspose(n1, vk);
1404  vk[1] = n3[1] * Trans.Weight();
1405  vk[2] = n3[2] * Trans.Weight();
1406  if (fe.GetMapType() == INTEGRAL)
1407  {
1408  double w = 1.0/Trans.Weight();
1409  for (int d = 0; d < 1; d++)
1410  {
1411  vk[d] *= w;
1412  }
1413  }
1414 
1415  for (int j = 0; j < shape.Size(); j++)
1416  {
1417  double s = shape(j);
1418  if (fabs(s) < 1e-12)
1419  {
1420  s = 0.0;
1421  }
1422  // Project scalar basis function multiplied by each coordinate
1423  // direction onto the transformed face normals
1424  for (int d = 0; d < vdim; d++)
1425  {
1426  I(k,j+d*shape.Size()) = s*vk[d];
1427  }
1428  }
1429  }
1430  }
1431  else
1432  {
1433  double vk[Geometry::MaxDim];
1434  DenseMatrix vshape(fe.GetDof(), fe.GetVDim());
1435 
1436  double * nk_ptr = const_cast<double*>(nk);
1437 
1438  I.SetSize(dof, fe.GetDof());
1439  for (int k = 0; k < dof; k++)
1440  {
1441  const IntegrationPoint &ip = Nodes.IntPoint(k);
1442 
1443  Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1444  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1445 
1446  Trans.SetIntPoint(&ip);
1447  // Transform RT face normals from reference to physical space
1448  // vk = adj(J)^T nk
1449  Trans.AdjugateJacobian().MultTranspose(n1, vk);
1450  // Compute fe basis functions in physical space
1451  fe.CalcVShape(Trans, vshape);
1452  // Project fe basis functions onto transformed face normals
1453  for (int j=0; j<vshape.Height(); j++)
1454  {
1455  I(k, j) = 0.0;
1456  I(k, j) += vshape(j, 0) * vk[0];
1457  if (vshape.Width() == 3)
1458  {
1459  I(k, j) += Trans.Weight() * vshape(j, 1) * n3(1);
1460  I(k, j) += Trans.Weight() * vshape(j, 2) * n3(2);
1461  }
1462  }
1463  }
1464  }
1465 }
1466 
1469  DenseMatrix &curl) const
1470 {
1471  DenseMatrix curl_shape(fe.GetDof(), fe.GetVDim());
1472  Vector curl_k(fe.GetDof());
1473 
1474  double * nk_ptr = const_cast<double*>(nk);
1475 
1476  curl.SetSize(dof, fe.GetDof());
1477  for (int k = 0; k < dof; k++)
1478  {
1479  fe.CalcCurlShape(Nodes.IntPoint(k), curl_shape);
1480  curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
1481  for (int j = 0; j < curl_k.Size(); j++)
1482  {
1483  curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1484  }
1485  }
1486 }
1487 
1488 const double RT_R2D_SegmentElement::nk[2] = { 0.,1.};
1489 
1491  const int ob_type)
1492  : VectorFiniteElement(1, Geometry::SEGMENT, p + 1, p + 1,
1493  H_DIV, FunctionSpace::Pk),
1494  dof2nk(dof),
1495  obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
1496 {
1497  // Override default dimension for VectorFiniteElements
1498  vdim = 2;
1499 
1500  const double *op = poly1d.OpenPoints(p, ob_type);
1501 
1502 #ifndef MFEM_THREAD_SAFE
1503  shape_ox.SetSize(p+1);
1504 #endif
1505 
1506  dof_map.SetSize(dof);
1507 
1508  int o = 0;
1509  // interior
1510  // z-components
1511  for (int i = 0; i <= p; i++)
1512  {
1513  Nodes.IntPoint(o).x = op[i];
1514  dof_map[i] = o; dof2nk[o++] = 0;
1515  }
1516 }
1517 
1519  DenseMatrix &shape) const
1520 {
1521  const int p = order;
1522 
1523 #ifdef MFEM_THREAD_SAFE
1524  Vector shape_ox(p);
1525 #endif
1526 
1527  obasis1d.Eval(ip.x, shape_ox);
1528 
1529  int o = 0;
1530  // z-components
1531  for (int i = 0; i <= p; i++)
1532  {
1533  int idx = dof_map[o++];
1534  shape(idx,0) = shape_ox(i);
1535  shape(idx,1) = 0.;
1536  }
1537 }
1538 
1540  DenseMatrix &shape) const
1541 {
1542  CalcVShape(Trans.GetIntPoint(), shape);
1543  const DenseMatrix & J = Trans.Jacobian();
1544  MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1545  "RT_R2D_SegmentElement cannot be embedded in "
1546  "2 or 3 dimensional spaces");
1547  for (int i=0; i<dof; i++)
1548  {
1549  shape(i, 0) *= J(0,0);
1550  }
1551  shape *= (1.0 / Trans.Weight());
1552 }
1553 
1555  Vector &div_shape) const
1556 {
1557  div_shape = 0.0;
1558 }
1559 
1560 void RT_R2D_SegmentElement::LocalInterpolation(const VectorFiniteElement &cfe,
1562  DenseMatrix &I) const
1563 {
1564  double vk[Geometry::MaxDim]; vk[1] = 0.0; vk[2] = 0.0;
1565  Vector xk(vk, dim);
1566  IntegrationPoint ip;
1567  DenseMatrix vshape(cfe.GetDof(), vdim);
1568 
1569  double * nk_ptr = const_cast<double*>(nk);
1570 
1571  I.SetSize(dof, vshape.Height());
1572 
1573  // assuming Trans is linear; this should be ok for all refinement types
1575  const DenseMatrix &adjJ = Trans.AdjugateJacobian();
1576  for (int k = 0; k < dof; k++)
1577  {
1578  Vector n2(&nk_ptr[dof2nk[k] * 2], 2);
1579 
1580  Trans.Transform(Nodes.IntPoint(k), xk);
1581  ip.Set3(vk);
1582  cfe.CalcVShape(ip, vshape);
1583  // xk = |J| J^{-t} n_k
1584  adjJ.MultTranspose(n2, vk);
1585  // I_k = vshape_k.adj(J)^t.n_k, k=1,...,dof
1586  for (int j = 0; j < vshape.Height(); j++)
1587  {
1588  double Ikj = 0.;
1589  /*
1590  for (int i = 0; i < dim; i++)
1591  {
1592  Ikj += vshape(j, i) * vk[i];
1593  }
1594  */
1595  Ikj += Trans.Weight() * vshape(j, 1) * n2(1);
1596  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1597  }
1598  }
1599 }
1600 
1602  const double *nk_fe)
1603  : VectorFiniteElement(2, G, Do, p + 1,
1604  H_DIV, FunctionSpace::Pk),
1605  nk(nk_fe),
1606  dof_map(dof),
1607  dof2nk(dof)
1608 {
1609  // Override default dimension for VectorFiniteElements
1610  vdim = 3;
1611 }
1612 
1614  DenseMatrix &shape) const
1615 {
1616  CalcVShape(Trans.GetIntPoint(), shape);
1617  const DenseMatrix & J = Trans.Jacobian();
1618  MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
1619  "RT_R2D_FiniteElement cannot be embedded in "
1620  "3 dimensional spaces");
1621  for (int i=0; i<dof; i++)
1622  {
1623  double sx = shape(i, 0);
1624  double sy = shape(i, 1);
1625  shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
1626  shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
1627  }
1628  shape *= (1.0 / Trans.Weight());
1629 }
1630 
1631 void
1632 RT_R2D_FiniteElement::LocalInterpolation(const VectorFiniteElement &cfe,
1633  ElementTransformation &Trans,
1634  DenseMatrix &I) const
1635 {
1636  double vk[Geometry::MaxDim]; vk[2] = 0.0;
1637  Vector xk(vk, dim);
1638  IntegrationPoint ip;
1639  DenseMatrix vshape(cfe.GetDof(), vdim);
1640 
1641  double * nk_ptr = const_cast<double*>(nk);
1642 
1643  I.SetSize(dof, vshape.Height());
1644 
1645  // assuming Trans is linear; this should be ok for all refinement types
1647  const DenseMatrix &adjJ = Trans.AdjugateJacobian();
1648  for (int k = 0; k < dof; k++)
1649  {
1650  Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
1651  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1652 
1653  Trans.Transform(Nodes.IntPoint(k), xk);
1654  ip.Set3(vk);
1655  cfe.CalcVShape(ip, vshape);
1656  // xk = |J| J^{-t} n_k
1657  adjJ.MultTranspose(n2, vk);
1658  // I_k = vshape_k.adj(J)^t.n_k, k=1,...,dof
1659  for (int j = 0; j < vshape.Height(); j++)
1660  {
1661  double Ikj = 0.;
1662  for (int i = 0; i < dim; i++)
1663  {
1664  Ikj += vshape(j, i) * vk[i];
1665  }
1666  Ikj += Trans.Weight() * vshape(j, 2) * n3(2);
1667  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1668  }
1669  }
1670 }
1671 
1673  DenseMatrix &R) const
1674 {
1675  double pt_data[Geometry::MaxDim];
1676  IntegrationPoint ip;
1677  Vector pt(pt_data, dim);
1678 
1679 #ifdef MFEM_THREAD_SAFE
1680  DenseMatrix vshape(dof, vdim);
1681 #endif
1682 
1683  double * nk_ptr = const_cast<double*>(nk);
1684 
1686  const DenseMatrix &J = Trans.Jacobian();
1687  const double weight = Trans.Weight();
1688  for (int j = 0; j < dof; j++)
1689  {
1690  Vector n2(&nk_ptr[dof2nk[j] * 3], 2);
1691  Vector n3(&nk_ptr[dof2nk[j] * 3], 3);
1692 
1693  InvertLinearTrans(Trans, Nodes.IntPoint(j), pt);
1694  ip.Set(pt_data, dim);
1695  if (Geometries.CheckPoint(geom_type, ip)) // do we need an epsilon here?
1696  {
1697  CalcVShape(ip, vshape);
1698  J.MultTranspose(n2, pt_data);
1699  pt /= weight;
1700  for (int k = 0; k < dof; k++)
1701  {
1702  double R_jk = 0.0;
1703  for (int d = 0; d < dim; d++)
1704  {
1705  R_jk += vshape(k,d)*pt_data[d];
1706  }
1707  R_jk += vshape(k,2) * n3(2);
1708  R(j,k) = R_jk;
1709  }
1710  }
1711  else
1712  {
1713  // Set the whole row to avoid valgrind warnings in R.Threshold().
1714  R.SetRow(j, infinity());
1715  }
1716  }
1717  R.Threshold(1e-12);
1718 }
1719 
1721  ElementTransformation &Trans,
1722  Vector &dofs) const
1723 {
1724  double data[3];
1725  Vector vk2(data, 2);
1726  Vector vk3(data, 3);
1727 
1728  double * nk_ptr = const_cast<double*>(nk);
1729 
1730  for (int k = 0; k < dof; k++)
1731  {
1732  Trans.SetIntPoint(&Nodes.IntPoint(k));
1733 
1734  vc.Eval(vk3, Trans, Nodes.IntPoint(k));
1735  // dof_k = nk^t adj(J) vk
1736  Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
1737  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1738 
1739  dofs(k) = Trans.AdjugateJacobian().InnerProduct(vk2, n2) +
1740  Trans.Weight() * vk3(2) * n3(2);
1741  }
1742 }
1743 
1745  ElementTransformation &Trans,
1746  DenseMatrix &I) const
1747 {
1748  if (fe.GetRangeType() == SCALAR)
1749  {
1750  double vk[Geometry::MaxDim];
1751  Vector shape(fe.GetDof());
1752 
1753  double * nk_ptr = const_cast<double*>(nk);
1754 
1755  I.SetSize(dof, vdim*fe.GetDof());
1756  for (int k = 0; k < dof; k++)
1757  {
1758  const IntegrationPoint &ip = Nodes.IntPoint(k);
1759 
1760  Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
1761  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1762 
1763  fe.CalcShape(ip, shape);
1764  Trans.SetIntPoint(&ip);
1765  // Transform RT face normals from reference to physical space
1766  // vk = adj(J)^T nk
1767  Trans.AdjugateJacobian().MultTranspose(n2, vk);
1768  vk[2] = n3[2] * Trans.Weight();
1769  if (fe.GetMapType() == INTEGRAL)
1770  {
1771  double w = 1.0/Trans.Weight();
1772  for (int d = 0; d < 2; d++)
1773  {
1774  vk[d] *= w;
1775  }
1776  }
1777 
1778  for (int j = 0; j < shape.Size(); j++)
1779  {
1780  double s = shape(j);
1781  if (fabs(s) < 1e-12)
1782  {
1783  s = 0.0;
1784  }
1785  // Project scalar basis function multiplied by each coordinate
1786  // direction onto the transformed face normals
1787  for (int d = 0; d < vdim; d++)
1788  {
1789  I(k,j+d*shape.Size()) = s*vk[d];
1790  }
1791  }
1792  }
1793  }
1794  else
1795  {
1796  double vk[Geometry::MaxDim];
1797  DenseMatrix vshape(fe.GetDof(), fe.GetVDim());
1798 
1799  double * nk_ptr = const_cast<double*>(nk);
1800 
1801  I.SetSize(dof, fe.GetDof());
1802  for (int k = 0; k < dof; k++)
1803  {
1804  const IntegrationPoint &ip = Nodes.IntPoint(k);
1805 
1806  Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
1807  Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1808 
1809  Trans.SetIntPoint(&ip);
1810  // Transform RT face normals from reference to physical space
1811  // vk = adj(J)^T nk
1812  Trans.AdjugateJacobian().MultTranspose(n2, vk);
1813  // Compute fe basis functions in physical space
1814  fe.CalcVShape(Trans, vshape);
1815  // Project fe basis functions onto transformed face normals
1816  for (int j=0; j<vshape.Height(); j++)
1817  {
1818  I(k, j) = 0.0;
1819  for (int i=0; i<2; i++)
1820  {
1821  I(k, j) += vshape(j, i) * vk[i];
1822  }
1823  if (vshape.Width() == 3)
1824  {
1825  I(k, j) += Trans.Weight() * vshape(j, 2) * n3(2);
1826  }
1827  }
1828  }
1829  }
1830 }
1831 
1833  ElementTransformation &Trans,
1834  DenseMatrix &curl) const
1835 {
1836  DenseMatrix curl_shape(fe.GetDof(), fe.GetVDim());
1837  Vector curl_k(fe.GetDof());
1838 
1839  double * nk_ptr = const_cast<double*>(nk);
1840 
1841  curl.SetSize(dof, fe.GetDof());
1842  for (int k = 0; k < dof; k++)
1843  {
1844  fe.CalcCurlShape(Nodes.IntPoint(k), curl_shape);
1845  curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
1846  for (int j = 0; j < curl_k.Size(); j++)
1847  {
1848  curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1849  }
1850  }
1851 }
1852 
1853 const double RT_R2D_TriangleElement::nk_t[12] =
1854 { 0.,-1.,0., 1.,1.,0., -1.,0.,0., 0.,0.,1. };
1855 
1857  : RT_R2D_FiniteElement(p, Geometry::TRIANGLE, ((p + 1)*(3 * p + 8))/2, nk_t),
1858  RT_FE(p),
1859  L2_FE(p)
1860 {
1861  L2_FE.SetMapType(INTEGRAL);
1862 
1863 #ifndef MFEM_THREAD_SAFE
1864  rt_shape.SetSize(RT_FE.GetDof(), 2);
1865  l2_shape.SetSize(L2_FE.GetDof());
1866  rt_dshape.SetSize(RT_FE.GetDof());
1867 #endif
1868 
1869  int o = 0;
1870  int r = 0;
1871  int l = 0;
1872 
1873  // Three edges
1874  for (int e=0; e<3; e++)
1875  {
1876  // Dofs in the plane
1877  for (int i=0; i<=p; i++)
1878  {
1879  dof_map[o] = r++; dof2nk[o++] = e;
1880  }
1881  }
1882 
1883  // Interior dofs in the plane
1884  for (int j = 0; j < p; j++)
1885  for (int i = 0; i + j < p; i++)
1886  {
1887  dof_map[o] = r++; dof2nk[o++] = 0;
1888  dof_map[o] = r++; dof2nk[o++] = 2;
1889  }
1890 
1891  // Interior z-directed dofs
1892  for (int j = 0; j <= p; j++)
1893  for (int i = 0; i + j <= p; i++)
1894  {
1895  dof_map[o] = -1 - l++; dof2nk[o++] = 3;
1896  }
1897 
1898  MFEM_VERIFY(r == RT_FE.GetDof(),
1899  "RT_R2D_Triangle incorrect number of RT dofs.");
1900  MFEM_VERIFY(l == L2_FE.GetDof(),
1901  "RT_R2D_Triangle incorrect number of L2 dofs.");
1902  MFEM_VERIFY(o == GetDof(),
1903  "RT_R2D_Triangle incorrect number of dofs.");
1904 
1905  const IntegrationRule & rt_Nodes = RT_FE.GetNodes();
1906  const IntegrationRule & l2_Nodes = L2_FE.GetNodes();
1907 
1908  for (int i=0; i<dof; i++)
1909  {
1910  int idx = dof_map[i];
1911  if (idx >= 0)
1912  {
1913  const IntegrationPoint & ip = rt_Nodes.IntPoint(idx);
1914  Nodes.IntPoint(i).Set2(ip.x, ip.y);
1915  }
1916  else
1917  {
1918  const IntegrationPoint & ip = l2_Nodes.IntPoint(-idx-1);
1919  Nodes.IntPoint(i).Set2(ip.x, ip.y);
1920  }
1921  }
1922 }
1923 
1925  DenseMatrix &shape) const
1926 {
1927 #ifdef MFEM_THREAD_SAFE
1928  DenseMatrix rt_shape(RT_FE.GetDof(), 2);
1929  Vector l2_shape(L2_FE.GetDof());
1930 #endif
1931 
1932  RT_FE.CalcVShape(ip, rt_shape);
1933  L2_FE.CalcShape(ip, l2_shape);
1934 
1935  for (int i=0; i<dof; i++)
1936  {
1937  int idx = dof_map[i];
1938  if (idx >= 0)
1939  {
1940  shape(i, 0) = rt_shape(idx, 0);
1941  shape(i, 1) = rt_shape(idx, 1);
1942  shape(i, 2) = 0.0;
1943  }
1944  else
1945  {
1946  shape(i, 0) = 0.0;
1947  shape(i, 1) = 0.0;
1948  shape(i, 2) = l2_shape(-idx-1);
1949  }
1950  }
1951 }
1952 
1954  Vector &div_shape) const
1955 {
1956 #ifdef MFEM_THREAD_SAFE
1957  Vector rt_dshape(RT_FE.GetDof());
1958 #endif
1959 
1960  RT_FE.CalcDivShape(ip, rt_dshape);
1961 
1962  for (int i=0; i<dof; i++)
1963  {
1964  int idx = dof_map[i];
1965  if (idx >= 0)
1966  {
1967  div_shape(i) = rt_dshape(idx);
1968  }
1969  else
1970  {
1971  div_shape(i) = 0.0;
1972  }
1973  }
1974 }
1975 
1976 const double RT_R2D_QuadrilateralElement::nk_q[15] =
1977 { 0., -1., 0., 1., 0., 0., 0., 1., 0., -1., 0., 0., 0., 0., 1. };
1978 
1980  const int cb_type,
1981  const int ob_type)
1982  : RT_R2D_FiniteElement(p, Geometry::SQUARE, (3*p + 5)*(p + 1), nk_q),
1983  cbasis1d(poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1984  obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
1985 {
1986  const double *cp = poly1d.ClosedPoints(p + 1, cb_type);
1987  const double *op = poly1d.OpenPoints(p, ob_type);
1988  const int dofx = (p + 1)*(p + 2);
1989  const int dofy = (p + 1)*(p + 2);
1990  const int dofxy = dofx + dofy;
1991 
1992 #ifndef MFEM_THREAD_SAFE
1993  shape_cx.SetSize(p + 2);
1994  shape_ox.SetSize(p + 1);
1995  shape_cy.SetSize(p + 2);
1996  shape_oy.SetSize(p + 1);
1997  dshape_cx.SetSize(p + 2);
1998  dshape_cy.SetSize(p + 2);
1999 #endif
2000 
2001  // edges
2002  int o = 0;
2003  for (int i = 0; i <= p; i++) // (0,1)
2004  {
2005  dof_map[dofx + i + 0*(p + 1)] = o++;
2006  }
2007  for (int i = 0; i <= p; i++) // (1,2)
2008  {
2009  dof_map[(p + 1) + i*(p + 2)] = o++;
2010  }
2011  for (int i = 0; i <= p; i++) // (2,3)
2012  {
2013  dof_map[dofx + (p - i) + (p + 1)*(p + 1)] = o++;
2014  }
2015  for (int i = 0; i <= p; i++) // (3,0)
2016  {
2017  dof_map[0 + (p - i)*(p + 2)] = o++;
2018  }
2019 
2020  // interior
2021  for (int j = 0; j <= p; j++) // x-components
2022  for (int i = 1; i <= p; i++)
2023  {
2024  dof_map[i + j*(p + 2)] = o++;
2025  }
2026  for (int j = 1; j <= p; j++) // y-components
2027  for (int i = 0; i <= p; i++)
2028  {
2029  dof_map[dofx + i + j*(p + 1)] = o++;
2030  }
2031  for (int j = 0; j <= p; j++) // z-components
2032  for (int i = 0; i <= p; i++)
2033  {
2034  dof_map[dofxy + i + j*(p + 1)] = o++;
2035  }
2036 
2037  // dof orientations
2038  // x-components
2039  for (int j = 0; j <= p; j++)
2040  for (int i = 0; i <= p/2; i++)
2041  {
2042  int idx = i + j*(p + 2);
2043  dof_map[idx] = -1 - dof_map[idx];
2044  }
2045  if (p%2 == 1)
2046  for (int j = p/2 + 1; j <= p; j++)
2047  {
2048  int idx = (p/2 + 1) + j*(p + 2);
2049  dof_map[idx] = -1 - dof_map[idx];
2050  }
2051  // y-components
2052  for (int j = 0; j <= p/2; j++)
2053  for (int i = 0; i <= p; i++)
2054  {
2055  int idx = dofx + i + j*(p + 1);
2056  dof_map[idx] = -1 - dof_map[idx];
2057  }
2058  if (p%2 == 1)
2059  for (int i = 0; i <= p/2; i++)
2060  {
2061  int idx = dofx + i + (p/2 + 1)*(p + 1);
2062  dof_map[idx] = -1 - dof_map[idx];
2063  }
2064 
2065  o = 0;
2066  for (int j = 0; j <= p; j++)
2067  for (int i = 0; i <= p + 1; i++)
2068  {
2069  int idx;
2070  if ((idx = dof_map[o++]) < 0)
2071  {
2072  idx = -1 - idx;
2073  dof2nk[idx] = 3;
2074  }
2075  else
2076  {
2077  dof2nk[idx] = 1;
2078  }
2079  Nodes.IntPoint(idx).Set2(cp[i], op[j]);
2080  }
2081  for (int j = 0; j <= p + 1; j++)
2082  for (int i = 0; i <= p; i++)
2083  {
2084  int idx;
2085  if ((idx = dof_map[o++]) < 0)
2086  {
2087  idx = -1 - idx;
2088  dof2nk[idx] = 0;
2089  }
2090  else
2091  {
2092  dof2nk[idx] = 2;
2093  }
2094  Nodes.IntPoint(idx).Set2(op[i], cp[j]);
2095  }
2096  for (int j = 0; j <= p; j++)
2097  for (int i = 0; i <= p; i++)
2098  {
2099  int idx = dof_map[o++];
2100  dof2nk[idx] = 4;
2101  Nodes.IntPoint(idx).Set2(op[i], op[j]);
2102  }
2103 }
2104 
2106  DenseMatrix &shape) const
2107 {
2108  const int pp1 = order;
2109 
2110 #ifdef MFEM_THREAD_SAFE
2111  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2112 #endif
2113 
2114  cbasis1d.Eval(ip.x, shape_cx);
2115  obasis1d.Eval(ip.x, shape_ox);
2116  cbasis1d.Eval(ip.y, shape_cy);
2117  obasis1d.Eval(ip.y, shape_oy);
2118 
2119  int o = 0;
2120  for (int j = 0; j < pp1; j++)
2121  for (int i = 0; i <= pp1; i++)
2122  {
2123  int idx, s;
2124  if ((idx = dof_map[o++]) < 0)
2125  {
2126  idx = -1 - idx, s = -1;
2127  }
2128  else
2129  {
2130  s = +1;
2131  }
2132  shape(idx,0) = s*shape_cx(i)*shape_oy(j);
2133  shape(idx,1) = 0.;
2134  shape(idx,2) = 0.;
2135  }
2136  for (int j = 0; j <= pp1; j++)
2137  for (int i = 0; i < pp1; i++)
2138  {
2139  int idx, s;
2140  if ((idx = dof_map[o++]) < 0)
2141  {
2142  idx = -1 - idx, s = -1;
2143  }
2144  else
2145  {
2146  s = +1;
2147  }
2148  shape(idx,0) = 0.;
2149  shape(idx,1) = s*shape_ox(i)*shape_cy(j);
2150  shape(idx,2) = 0.;
2151  }
2152  for (int j = 0; j < pp1; j++)
2153  for (int i = 0; i < pp1; i++)
2154  {
2155  int idx = dof_map[o++];
2156  shape(idx,0) = 0.;
2157  shape(idx,1) = 0.;
2158  shape(idx,2) = shape_ox(i)*shape_oy(j);
2159  }
2160 }
2161 
2163  Vector &divshape) const
2164 {
2165  const int pp1 = order;
2166 
2167 #ifdef MFEM_THREAD_SAFE
2168  Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2169  Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
2170 #endif
2171 
2172  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
2173  obasis1d.Eval(ip.x, shape_ox);
2174  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
2175  obasis1d.Eval(ip.y, shape_oy);
2176 
2177  int o = 0;
2178  for (int j = 0; j < pp1; j++)
2179  for (int i = 0; i <= pp1; i++)
2180  {
2181  int idx, s;
2182  if ((idx = dof_map[o++]) < 0)
2183  {
2184  idx = -1 - idx, s = -1;
2185  }
2186  else
2187  {
2188  s = +1;
2189  }
2190  divshape(idx) = s*dshape_cx(i)*shape_oy(j);
2191  }
2192  for (int j = 0; j <= pp1; j++)
2193  for (int i = 0; i < pp1; i++)
2194  {
2195  int idx, s;
2196  if ((idx = dof_map[o++]) < 0)
2197  {
2198  idx = -1 - idx, s = -1;
2199  }
2200  else
2201  {
2202  s = +1;
2203  }
2204  divshape(idx) = s*shape_ox(i)*dshape_cy(j);
2205  }
2206  for (int j = 0; j < pp1; j++)
2207  for (int i = 0; i < pp1; i++)
2208  {
2209  int idx = dof_map[o++];
2210  divshape(idx) = 0.;
2211  }
2212 }
2213 
2214 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
void Set(const double *p, const int dim)
Definition: intrules.hpp:37
const DenseMatrix & AdjugateJacobian()
Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoin...
Definition: eltrans.hpp:135
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:946
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_base.cpp:39
void ScaleIntegrated(bool scale_integrated_)
Set whether the &quot;integrated&quot; basis should be scaled by the subcell sizes. Has no effect for non-integ...
Definition: fe_base.cpp:1879
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
bool IsIntegratedType() const
Returns true if the basis is &quot;integrated&quot;, false otherwise.
Definition: fe_base.hpp:1022
Base class for vector Coefficients that optionally depend on time and space.
void SetRow(int r, const double *row)
Definition: densemat.cpp:2139
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:1268
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_rt.cpp:1467
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:200
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:1145
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:2105
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
int dim
Dimension of reference space.
Definition: fe_base.hpp:238
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:297
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_rt.cpp:1832
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:802
DenseMatrix vshape
Definition: fe_base.hpp:250
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe_base.hpp:349
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition: fe_l2.cpp:39
void Factor()
Factor the current DenseMatrix, *a.
Definition: densemat.cpp:3878
Geometry::Type geom_type
Geometry::Type of the reference element.
Definition: fe_base.hpp:241
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
Intermediate class for finite elements whose basis functions return vector values.
Definition: fe_base.hpp:786
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:1924
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:71
const double * nk
Definition: fe_rt.hpp:423
RT_WedgeElement(const int p)
Definition: fe_rt.cpp:1027
const double * ClosedPoints(const int p, const int btype=BasisType::GaussLobatto)
Get coordinates of a closed (GaussLegendre) set of points if degree p.
Definition: fe_base.hpp:1070
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe_rt.cpp:635
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:982
int vdim
Vector dimension of vector-valued basis functions.
Definition: fe_base.hpp:239
static void CalcBasis(const int p, const double x, double *u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
Definition: fe_base.hpp:1086
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:1322
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
Array< int > dof_map
Definition: fe_rt.hpp:424
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:201
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.hpp:668
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in physical space at the point ...
Definition: fe_rt.cpp:1613
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
RT_R1D_SegmentElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R1D_SegmentElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_rt.cpp:1214
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &div_shape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:1554
RT_R2D_FiniteElement(int p, Geometry::Type G, int Do, const double *nk_fe)
Definition: fe_rt.cpp:1601
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:2162
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the &quot;integrated&quot; basis type using pre-computed closed basis derivatives. ...
Definition: fe_base.cpp:1854
void Set2(const double x1, const double x2)
Definition: intrules.hpp:80
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:468
int GetVDim()
Returns dimension of the vector.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition: fe_l2.cpp:615
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:141
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:556
const double * OpenPoints(const int p, const int btype=BasisType::GaussLegendre)
Get coordinates of an open (GaussLegendre) set of points if degree p.
Definition: fe_base.hpp:1065
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:1180
static const int MaxDim
Definition: geom.hpp:43
void Threshold(double eps)
Replace small entries, abs(a_ij) &lt;= eps, with zero.
Definition: densemat.cpp:2169
RT_R2D_TriangleElement(const int p)
Construct the RT_R2D_TriangleElement of order p.
Definition: fe_rt.cpp:1856
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:70
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:1518
IntegrationRule Nodes
Definition: fe_base.hpp:248
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe_rt.cpp:255
RT_R2D_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_rt.cpp:1979
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
RT_R2D_SegmentElement(const int p, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R2D_SegmentElement of order p and open BasisType ob_type.
Definition: fe_rt.cpp:1490
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:3924
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
Definition: fe_base.cpp:64
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Definition: fe_rt.cpp:1355
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe_rt.cpp:1672
int dof
Number of degrees of freedom.
Definition: fe_base.hpp:245
RT_TetrahedronElement(const int p)
Construct the RT_TetrahedronElement of order p.
Definition: fe_rt.cpp:844
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition: fe_h1.cpp:59
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Definition: fe_rt.cpp:1720
RT_TriangleElement(const int p)
Construct the RT_TriangleElement of order p.
Definition: fe_rt.cpp:691
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition: fe_h1.cpp:40
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
Definition: fe_base.cpp:1622
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
Definition: fe_base.cpp:611
Vector data type.
Definition: vector.hpp:60
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_rt.cpp:769
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition: geom.cpp:405
Describes the function space on each element.
Definition: fe_base.hpp:218
RefCoord s[3]
RT_HexahedronElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_HexahedronElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_rt.cpp:302
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
RT_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_rt.cpp:25
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition: fe_base.hpp:340
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
int order
Order/degree of the shape functions.
Definition: fe_base.hpp:245
Poly_1D poly1d
Definition: fe.cpp:28
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_rt.cpp:1953