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