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