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