MFEM  v4.5.2
Finite element discretization library
fe_nd.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 // Nedelec Finite Element classes
13 
14 #include "fe_nd.hpp"
15 #include "../coefficient.hpp"
16 
17 namespace mfem
18 {
19 
20 using namespace std;
21 
22 const double ND_HexahedronElement::tk[18] =
23 { 1.,0.,0., 0.,1.,0., 0.,0.,1., -1.,0.,0., 0.,-1.,0., 0.,0.,-1. };
24 
26  const int cb_type, const int ob_type)
27  : VectorTensorFiniteElement(3, 3*p*(p + 1)*(p + 1), p, cb_type, ob_type,
28  H_CURL, DofMapType::L2_DOF_MAP),
29  dof2tk(dof), cp(poly1d.ClosedPoints(p, cb_type))
30 {
31  if (obasis1d.IsIntegratedType()) { is_nodal = false; }
32 
34 
35  const double *op = poly1d.OpenPoints(p - 1, ob_type);
36  const int dof3 = dof/3;
37 
38 #ifndef MFEM_THREAD_SAFE
39  shape_cx.SetSize(p + 1);
40  shape_ox.SetSize(p);
41  shape_cy.SetSize(p + 1);
42  shape_oy.SetSize(p);
43  shape_cz.SetSize(p + 1);
44  shape_oz.SetSize(p);
45  dshape_cx.SetSize(p + 1);
46  dshape_cy.SetSize(p + 1);
47  dshape_cz.SetSize(p + 1);
48 #endif
49 
50  // edges
51  int o = 0;
52  for (int i = 0; i < p; i++) // (0,1)
53  {
54  dof_map[0*dof3 + i + (0 + 0*(p + 1))*p] = o++;
55  }
56  for (int i = 0; i < p; i++) // (1,2)
57  {
58  dof_map[1*dof3 + p + (i + 0*p)*(p + 1)] = o++;
59  }
60  for (int i = 0; i < p; i++) // (3,2)
61  {
62  dof_map[0*dof3 + i + (p + 0*(p + 1))*p] = o++;
63  }
64  for (int i = 0; i < p; i++) // (0,3)
65  {
66  dof_map[1*dof3 + 0 + (i + 0*p)*(p + 1)] = o++;
67  }
68  for (int i = 0; i < p; i++) // (4,5)
69  {
70  dof_map[0*dof3 + i + (0 + p*(p + 1))*p] = o++;
71  }
72  for (int i = 0; i < p; i++) // (5,6)
73  {
74  dof_map[1*dof3 + p + (i + p*p)*(p + 1)] = o++;
75  }
76  for (int i = 0; i < p; i++) // (7,6)
77  {
78  dof_map[0*dof3 + i + (p + p*(p + 1))*p] = o++;
79  }
80  for (int i = 0; i < p; i++) // (4,7)
81  {
82  dof_map[1*dof3 + 0 + (i + p*p)*(p + 1)] = o++;
83  }
84  for (int i = 0; i < p; i++) // (0,4)
85  {
86  dof_map[2*dof3 + 0 + (0 + i*(p + 1))*(p + 1)] = o++;
87  }
88  for (int i = 0; i < p; i++) // (1,5)
89  {
90  dof_map[2*dof3 + p + (0 + i*(p + 1))*(p + 1)] = o++;
91  }
92  for (int i = 0; i < p; i++) // (2,6)
93  {
94  dof_map[2*dof3 + p + (p + i*(p + 1))*(p + 1)] = o++;
95  }
96  for (int i = 0; i < p; i++) // (3,7)
97  {
98  dof_map[2*dof3 + 0 + (p + i*(p + 1))*(p + 1)] = o++;
99  }
100 
101  // faces
102  // (3,2,1,0) -- bottom
103  for (int j = 1; j < p; j++) // x - components
104  for (int i = 0; i < p; i++)
105  {
106  dof_map[0*dof3 + i + ((p - j) + 0*(p + 1))*p] = o++;
107  }
108  for (int j = 0; j < p; j++) // y - components
109  for (int i = 1; i < p; i++)
110  {
111  dof_map[1*dof3 + i + ((p - 1 - j) + 0*p)*(p + 1)] = -1 - (o++);
112  }
113  // (0,1,5,4) -- front
114  for (int k = 1; k < p; k++) // x - components
115  for (int i = 0; i < p; i++)
116  {
117  dof_map[0*dof3 + i + (0 + k*(p + 1))*p] = o++;
118  }
119  for (int k = 0; k < p; k++) // z - components
120  for (int i = 1; i < p; i++ )
121  {
122  dof_map[2*dof3 + i + (0 + k*(p + 1))*(p + 1)] = o++;
123  }
124  // (1,2,6,5) -- right
125  for (int k = 1; k < p; k++) // y - components
126  for (int j = 0; j < p; j++)
127  {
128  dof_map[1*dof3 + p + (j + k*p)*(p + 1)] = o++;
129  }
130  for (int k = 0; k < p; k++) // z - components
131  for (int j = 1; j < p; j++)
132  {
133  dof_map[2*dof3 + p + (j + k*(p + 1))*(p + 1)] = o++;
134  }
135  // (2,3,7,6) -- back
136  for (int k = 1; k < p; k++) // x - components
137  for (int i = 0; i < p; i++)
138  {
139  dof_map[0*dof3 + (p - 1 - i) + (p + k*(p + 1))*p] = -1 - (o++);
140  }
141  for (int k = 0; k < p; k++) // z - components
142  for (int i = 1; i < p; i++)
143  {
144  dof_map[2*dof3 + (p - i) + (p + k*(p + 1))*(p + 1)] = o++;
145  }
146  // (3,0,4,7) -- left
147  for (int k = 1; k < p; k++) // y - components
148  for (int j = 0; j < p; j++)
149  {
150  dof_map[1*dof3 + 0 + ((p - 1 - j) + k*p)*(p + 1)] = -1 - (o++);
151  }
152  for (int k = 0; k < p; k++) // z - components
153  for (int j = 1; j < p; j++)
154  {
155  dof_map[2*dof3 + 0 + ((p - j) + k*(p + 1))*(p + 1)] = o++;
156  }
157  // (4,5,6,7) -- top
158  for (int j = 1; j < p; j++) // x - components
159  for (int i = 0; i < p; i++)
160  {
161  dof_map[0*dof3 + i + (j + p*(p + 1))*p] = o++;
162  }
163  for (int j = 0; j < p; j++) // y - components
164  for (int i = 1; i < p; i++)
165  {
166  dof_map[1*dof3 + i + (j + p*p)*(p + 1)] = o++;
167  }
168 
169  // interior
170  // x-components
171  for (int k = 1; k < p; k++)
172  for (int j = 1; j < p; j++)
173  for (int i = 0; i < p; i++)
174  {
175  dof_map[0*dof3 + i + (j + k*(p + 1))*p] = o++;
176  }
177  // y-components
178  for (int k = 1; k < p; k++)
179  for (int j = 0; j < p; j++)
180  for (int i = 1; i < p; i++)
181  {
182  dof_map[1*dof3 + i + (j + k*p)*(p + 1)] = o++;
183  }
184  // z-components
185  for (int k = 0; k < p; k++)
186  for (int j = 1; j < p; j++)
187  for (int i = 1; i < p; i++)
188  {
189  dof_map[2*dof3 + i + (j + k*(p + 1))*(p + 1)] = o++;
190  }
191 
192  // set dof2tk and Nodes
193  o = 0;
194  // x-components
195  for (int k = 0; k <= p; k++)
196  for (int j = 0; j <= p; j++)
197  for (int i = 0; i < p; i++)
198  {
199  int idx;
200  if ((idx = dof_map[o++]) < 0)
201  {
202  dof2tk[idx = -1 - idx] = 3;
203  }
204  else
205  {
206  dof2tk[idx] = 0;
207  }
208  Nodes.IntPoint(idx).Set3(op[i], cp[j], cp[k]);
209  }
210  // y-components
211  for (int k = 0; k <= p; k++)
212  for (int j = 0; j < p; j++)
213  for (int i = 0; i <= p; i++)
214  {
215  int idx;
216  if ((idx = dof_map[o++]) < 0)
217  {
218  dof2tk[idx = -1 - idx] = 4;
219  }
220  else
221  {
222  dof2tk[idx] = 1;
223  }
224  Nodes.IntPoint(idx).Set3(cp[i], op[j], cp[k]);
225  }
226  // z-components
227  for (int k = 0; k < p; k++)
228  for (int j = 0; j <= p; j++)
229  for (int i = 0; i <= p; i++)
230  {
231  int idx;
232  if ((idx = dof_map[o++]) < 0)
233  {
234  dof2tk[idx = -1 - idx] = 5;
235  }
236  else
237  {
238  dof2tk[idx] = 2;
239  }
240  Nodes.IntPoint(idx).Set3(cp[i], cp[j], op[k]);
241  }
242 }
243 
246  Vector &dofs) const
247 {
248  MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
249  double vk[Geometry::MaxDim];
250  Vector xk(vk, vc.GetVDim());
251 
253  const int nqpt = ir.GetNPoints();
254 
255  IntegrationPoint ip3d;
256 
257  int o = 0;
258  for (int c = 0; c < 3; ++c) // loop over x, y, z components
259  {
260  const int im = c == 0 ? order - 1 : order;
261  const int jm = c == 1 ? order - 1 : order;
262  const int km = c == 2 ? order - 1 : order;
263 
264  for (int k = 0; k <= km; k++)
265  for (int j = 0; j <= jm; j++)
266  for (int i = 0; i <= im; i++)
267  {
268  int idx;
269  if ((idx = dof_map[o++]) < 0)
270  {
271  idx = -1 - idx;
272  }
273 
274  const int id1 = c == 0 ? i : (c == 1 ? j : k);
275  const double h = cp[id1+1] - cp[id1];
276 
277  double val = 0.0;
278 
279  for (int q = 0; q < nqpt; q++)
280  {
281  const IntegrationPoint &ip1d = ir.IntPoint(q);
282 
283  if (c == 0)
284  {
285  ip3d.Set3(cp[i] + (h*ip1d.x), cp[j], cp[k]);
286  }
287  else if (c == 1)
288  {
289  ip3d.Set3(cp[i], cp[j] + (h*ip1d.x), cp[k]);
290  }
291  else
292  {
293  ip3d.Set3(cp[i], cp[j], cp[k] + (h*ip1d.x));
294  }
295 
296  Trans.SetIntPoint(&ip3d);
297  vc.Eval(xk, Trans, ip3d);
298 
299  // xk^t J tk
300  const double ipval = Trans.Jacobian().InnerProduct(tk + dof2tk[idx]*dim, vk);
301  val += ip1d.weight * ipval;
302  }
303 
304  dofs(idx) = val*h;
305  }
306  }
307 }
308 
310  DenseMatrix &shape) const
311 {
312  const int p = order;
313 
314 #ifdef MFEM_THREAD_SAFE
315  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
316  Vector shape_cz(p + 1), shape_oz(p);
317 #endif
318 
320  {
321 #ifdef MFEM_THREAD_SAFE
322  Vector dshape_cx(p + 1), dshape_cy(p + 1), dshape_cz(p + 1);
323 #endif
324  basis1d.Eval(ip.x, shape_cx, dshape_cx);
325  basis1d.Eval(ip.y, shape_cy, dshape_cy);
326  basis1d.Eval(ip.z, shape_cz, dshape_cz);
327  obasis1d.ScaleIntegrated(false);
328  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
329  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
330  obasis1d.EvalIntegrated(dshape_cz, shape_oz);
331  }
332  else
333  {
334  basis1d.Eval(ip.x, shape_cx);
335  basis1d.Eval(ip.y, shape_cy);
336  basis1d.Eval(ip.z, shape_cz);
337  obasis1d.Eval(ip.x, shape_ox);
338  obasis1d.Eval(ip.y, shape_oy);
339  obasis1d.Eval(ip.z, shape_oz);
340  }
341 
342  int o = 0;
343  // x-components
344  for (int k = 0; k <= p; k++)
345  for (int j = 0; j <= p; j++)
346  for (int i = 0; i < p; i++)
347  {
348  int idx, s;
349  if ((idx = dof_map[o++]) < 0)
350  {
351  idx = -1 - idx, s = -1;
352  }
353  else
354  {
355  s = +1;
356  }
357  shape(idx,0) = s*shape_ox(i)*shape_cy(j)*shape_cz(k);
358  shape(idx,1) = 0.;
359  shape(idx,2) = 0.;
360  }
361  // y-components
362  for (int k = 0; k <= p; k++)
363  for (int j = 0; j < p; j++)
364  for (int i = 0; i <= p; i++)
365  {
366  int idx, s;
367  if ((idx = dof_map[o++]) < 0)
368  {
369  idx = -1 - idx, s = -1;
370  }
371  else
372  {
373  s = +1;
374  }
375  shape(idx,0) = 0.;
376  shape(idx,1) = s*shape_cx(i)*shape_oy(j)*shape_cz(k);
377  shape(idx,2) = 0.;
378  }
379  // z-components
380  for (int k = 0; k < p; k++)
381  for (int j = 0; j <= p; j++)
382  for (int i = 0; i <= p; i++)
383  {
384  int idx, s;
385  if ((idx = dof_map[o++]) < 0)
386  {
387  idx = -1 - idx, s = -1;
388  }
389  else
390  {
391  s = +1;
392  }
393  shape(idx,0) = 0.;
394  shape(idx,1) = 0.;
395  shape(idx,2) = s*shape_cx(i)*shape_cy(j)*shape_oz(k);
396  }
397 }
398 
400  DenseMatrix &curl_shape) const
401 {
402  const int p = order;
403 
404 #ifdef MFEM_THREAD_SAFE
405  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
406  Vector shape_cz(p + 1), shape_oz(p);
407  Vector dshape_cx(p + 1), dshape_cy(p + 1), dshape_cz(p + 1);
408 #endif
409 
410  basis1d.Eval(ip.x, shape_cx, dshape_cx);
411  basis1d.Eval(ip.y, shape_cy, dshape_cy);
412  basis1d.Eval(ip.z, shape_cz, dshape_cz);
414  {
415  obasis1d.ScaleIntegrated(false);
416  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
417  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
418  obasis1d.EvalIntegrated(dshape_cz, shape_oz);
419  }
420  else
421  {
422  obasis1d.Eval(ip.x, shape_ox);
423  obasis1d.Eval(ip.y, shape_oy);
424  obasis1d.Eval(ip.z, shape_oz);
425  }
426 
427  int o = 0;
428  // x-components
429  for (int k = 0; k <= p; k++)
430  for (int j = 0; j <= p; j++)
431  for (int i = 0; i < p; i++)
432  {
433  int idx, s;
434  if ((idx = dof_map[o++]) < 0)
435  {
436  idx = -1 - idx, s = -1;
437  }
438  else
439  {
440  s = +1;
441  }
442  curl_shape(idx,0) = 0.;
443  curl_shape(idx,1) = s*shape_ox(i)* shape_cy(j)*dshape_cz(k);
444  curl_shape(idx,2) = -s*shape_ox(i)*dshape_cy(j)* shape_cz(k);
445  }
446  // y-components
447  for (int k = 0; k <= p; k++)
448  for (int j = 0; j < p; j++)
449  for (int i = 0; i <= p; i++)
450  {
451  int idx, s;
452  if ((idx = dof_map[o++]) < 0)
453  {
454  idx = -1 - idx, s = -1;
455  }
456  else
457  {
458  s = +1;
459  }
460  curl_shape(idx,0) = -s* shape_cx(i)*shape_oy(j)*dshape_cz(k);
461  curl_shape(idx,1) = 0.;
462  curl_shape(idx,2) = s*dshape_cx(i)*shape_oy(j)* shape_cz(k);
463  }
464  // z-components
465  for (int k = 0; k < p; k++)
466  for (int j = 0; j <= p; j++)
467  for (int i = 0; i <= p; i++)
468  {
469  int idx, s;
470  if ((idx = dof_map[o++]) < 0)
471  {
472  idx = -1 - idx, s = -1;
473  }
474  else
475  {
476  s = +1;
477  }
478  curl_shape(idx,0) = s* shape_cx(i)*dshape_cy(j)*shape_oz(k);
479  curl_shape(idx,1) = -s*dshape_cx(i)* shape_cy(j)*shape_oz(k);
480  curl_shape(idx,2) = 0.;
481  }
482 }
483 
484 const double ND_QuadrilateralElement::tk[8] =
485 { 1.,0., 0.,1., -1.,0., 0.,-1. };
486 
488  const int cb_type,
489  const int ob_type)
490  : VectorTensorFiniteElement(2, 2*p*(p + 1), p, cb_type, ob_type,
491  H_CURL, DofMapType::L2_DOF_MAP),
492  dof2tk(dof),
493  cp(poly1d.ClosedPoints(p, cb_type))
494 {
495  if (obasis1d.IsIntegratedType()) { is_nodal = false; }
496 
498 
499  const double *op = poly1d.OpenPoints(p - 1, ob_type);
500  const int dof2 = dof/2;
501 
502 #ifndef MFEM_THREAD_SAFE
503  shape_cx.SetSize(p + 1);
504  shape_ox.SetSize(p);
505  shape_cy.SetSize(p + 1);
506  shape_oy.SetSize(p);
507  dshape_cx.SetSize(p + 1);
508  dshape_cy.SetSize(p + 1);
509 #endif
510 
511  // edges
512  int o = 0;
513  for (int i = 0; i < p; i++) // (0,1)
514  {
515  dof_map[0*dof2 + i + 0*p] = o++;
516  }
517  for (int j = 0; j < p; j++) // (1,2)
518  {
519  dof_map[1*dof2 + p + j*(p + 1)] = o++;
520  }
521  for (int i = 0; i < p; i++) // (2,3)
522  {
523  dof_map[0*dof2 + (p - 1 - i) + p*p] = -1 - (o++);
524  }
525  for (int j = 0; j < p; j++) // (3,0)
526  {
527  dof_map[1*dof2 + 0 + (p - 1 - j)*(p + 1)] = -1 - (o++);
528  }
529 
530  // interior
531  // x-components
532  for (int j = 1; j < p; j++)
533  for (int i = 0; i < p; i++)
534  {
535  dof_map[0*dof2 + i + j*p] = o++;
536  }
537  // y-components
538  for (int j = 0; j < p; j++)
539  for (int i = 1; i < p; i++)
540  {
541  dof_map[1*dof2 + i + j*(p + 1)] = o++;
542  }
543 
544  // set dof2tk and Nodes
545  o = 0;
546  // x-components
547  for (int j = 0; j <= p; j++)
548  for (int i = 0; i < p; i++)
549  {
550  int idx;
551  if ((idx = dof_map[o++]) < 0)
552  {
553  dof2tk[idx = -1 - idx] = 2;
554  }
555  else
556  {
557  dof2tk[idx] = 0;
558  }
559  Nodes.IntPoint(idx).Set2(op[i], cp[j]);
560  }
561  // y-components
562  for (int j = 0; j < p; j++)
563  for (int i = 0; i <= p; i++)
564  {
565  int idx;
566  if ((idx = dof_map[o++]) < 0)
567  {
568  dof2tk[idx = -1 - idx] = 3;
569  }
570  else
571  {
572  dof2tk[idx] = 1;
573  }
574  Nodes.IntPoint(idx).Set2(cp[i], op[j]);
575  }
576 }
577 
580  Vector &dofs) const
581 {
582  MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
583  double vk[Geometry::MaxDim];
584  Vector xk(vk, vc.GetVDim());
585 
587  const int nqpt = ir.GetNPoints();
588 
589  IntegrationPoint ip2d;
590 
591  int o = 0;
592  // x-components
593  for (int j = 0; j <= order; j++)
594  for (int i = 0; i < order; i++)
595  {
596  int idx;
597  if ((idx = dof_map[o++]) < 0)
598  {
599  idx = -1 - idx;
600  }
601 
602  const double h = cp[i+1] - cp[i];
603 
604  double val = 0.0;
605 
606  for (int k = 0; k < nqpt; k++)
607  {
608  const IntegrationPoint &ip1d = ir.IntPoint(k);
609 
610  ip2d.Set2(cp[i] + (h*ip1d.x), cp[j]);
611 
612  Trans.SetIntPoint(&ip2d);
613  vc.Eval(xk, Trans, ip2d);
614 
615  // xk^t J tk
616  const double ipval = Trans.Jacobian().InnerProduct(tk + dof2tk[idx]*dim, vk);
617  val += ip1d.weight * ipval;
618  }
619 
620  dofs(idx) = val*h;
621  }
622  // y-components
623  for (int j = 0; j < order; j++)
624  for (int i = 0; i <= order; i++)
625  {
626  int idx;
627  if ((idx = dof_map[o++]) < 0)
628  {
629  idx = -1 - idx;
630  }
631 
632  const double h = cp[j+1] - cp[j];
633 
634  double val = 0.0;
635 
636  for (int k = 0; k < nqpt; k++)
637  {
638  const IntegrationPoint &ip1d = ir.IntPoint(k);
639 
640  ip2d.Set2(cp[i], cp[j] + (h*ip1d.x));
641 
642  Trans.SetIntPoint(&ip2d);
643  vc.Eval(xk, Trans, ip2d);
644 
645  // xk^t J tk
646  const double ipval = Trans.Jacobian().InnerProduct(tk + dof2tk[idx]*dim, vk);
647  val += ip1d.weight * ipval;
648  }
649 
650  dofs(idx) = val*h;
651  }
652 }
653 
655  DenseMatrix &shape) const
656 {
657  const int p = order;
658 
659 #ifdef MFEM_THREAD_SAFE
660  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
661 #endif
662 
664  {
665 #ifdef MFEM_THREAD_SAFE
666  Vector dshape_cx(p + 1), dshape_cy(p + 1);
667 #endif
668  basis1d.Eval(ip.x, shape_cx, dshape_cx);
669  basis1d.Eval(ip.y, shape_cy, dshape_cy);
670  obasis1d.ScaleIntegrated(false);
671  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
672  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
673  }
674  else
675  {
676  basis1d.Eval(ip.x, shape_cx);
677  basis1d.Eval(ip.y, shape_cy);
678  obasis1d.Eval(ip.x, shape_ox);
679  obasis1d.Eval(ip.y, shape_oy);
680  }
681 
682  int o = 0;
683  // x-components
684  for (int j = 0; j <= p; j++)
685  for (int i = 0; i < p; i++)
686  {
687  int idx, s;
688  if ((idx = dof_map[o++]) < 0)
689  {
690  idx = -1 - idx, s = -1;
691  }
692  else
693  {
694  s = +1;
695  }
696  shape(idx,0) = s*shape_ox(i)*shape_cy(j);
697  shape(idx,1) = 0.;
698  }
699  // y-components
700  for (int j = 0; j < p; j++)
701  for (int i = 0; i <= p; i++)
702  {
703  int idx, s;
704  if ((idx = dof_map[o++]) < 0)
705  {
706  idx = -1 - idx, s = -1;
707  }
708  else
709  {
710  s = +1;
711  }
712  shape(idx,0) = 0.;
713  shape(idx,1) = s*shape_cx(i)*shape_oy(j);
714  }
715 }
716 
718  DenseMatrix &curl_shape) const
719 {
720  const int p = order;
721 
722 #ifdef MFEM_THREAD_SAFE
723  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
724  Vector dshape_cx(p + 1), dshape_cy(p + 1);
725 #endif
726 
727  basis1d.Eval(ip.x, shape_cx, dshape_cx);
728  basis1d.Eval(ip.y, shape_cy, dshape_cy);
730  {
731  obasis1d.ScaleIntegrated(false);
732  obasis1d.EvalIntegrated(dshape_cx, shape_ox);
733  obasis1d.EvalIntegrated(dshape_cy, shape_oy);
734  }
735  else
736  {
737  obasis1d.Eval(ip.x, shape_ox);
738  obasis1d.Eval(ip.y, shape_oy);
739  }
740 
741  int o = 0;
742  // x-components
743  for (int j = 0; j <= p; j++)
744  for (int i = 0; i < p; i++)
745  {
746  int idx, s;
747  if ((idx = dof_map[o++]) < 0)
748  {
749  idx = -1 - idx, s = -1;
750  }
751  else
752  {
753  s = +1;
754  }
755  curl_shape(idx,0) = -s*shape_ox(i)*dshape_cy(j);
756  }
757  // y-components
758  for (int j = 0; j < p; j++)
759  for (int i = 0; i <= p; i++)
760  {
761  int idx, s;
762  if ((idx = dof_map[o++]) < 0)
763  {
764  idx = -1 - idx, s = -1;
765  }
766  else
767  {
768  s = +1;
769  }
770  curl_shape(idx,0) = s*dshape_cx(i)*shape_oy(j);
771  }
772 }
773 
774 
775 const double ND_TetrahedronElement::tk[18] =
776 { 1.,0.,0., 0.,1.,0., 0.,0.,1., -1.,1.,0., -1.,0.,1., 0.,-1.,1. };
777 
778 const double ND_TetrahedronElement::c = 1./4.;
779 
781  : VectorFiniteElement(3, Geometry::TETRAHEDRON, p*(p + 2)*(p + 3)/2, p,
782  H_CURL, FunctionSpace::Pk), dof2tk(dof)
783 {
784  const double *eop = poly1d.OpenPoints(p - 1);
785  const double *fop = (p > 1) ? poly1d.OpenPoints(p - 2) : NULL;
786  const double *iop = (p > 2) ? poly1d.OpenPoints(p - 3) : NULL;
787 
788  const int pm1 = p - 1, pm2 = p - 2, pm3 = p - 3;
789 
790 #ifndef MFEM_THREAD_SAFE
791  shape_x.SetSize(p);
792  shape_y.SetSize(p);
793  shape_z.SetSize(p);
794  shape_l.SetSize(p);
795  dshape_x.SetSize(p);
796  dshape_y.SetSize(p);
797  dshape_z.SetSize(p);
798  dshape_l.SetSize(p);
799  u.SetSize(dof, dim);
800 #else
801  Vector shape_x(p), shape_y(p), shape_z(p), shape_l(p);
802 #endif
803 
804  int o = 0;
805  // edges
806  for (int i = 0; i < p; i++) // (0,1)
807  {
808  Nodes.IntPoint(o).Set3(eop[i], 0., 0.);
809  dof2tk[o++] = 0;
810  }
811  for (int i = 0; i < p; i++) // (0,2)
812  {
813  Nodes.IntPoint(o).Set3(0., eop[i], 0.);
814  dof2tk[o++] = 1;
815  }
816  for (int i = 0; i < p; i++) // (0,3)
817  {
818  Nodes.IntPoint(o).Set3(0., 0., eop[i]);
819  dof2tk[o++] = 2;
820  }
821  for (int i = 0; i < p; i++) // (1,2)
822  {
823  Nodes.IntPoint(o).Set3(eop[pm1-i], eop[i], 0.);
824  dof2tk[o++] = 3;
825  }
826  for (int i = 0; i < p; i++) // (1,3)
827  {
828  Nodes.IntPoint(o).Set3(eop[pm1-i], 0., eop[i]);
829  dof2tk[o++] = 4;
830  }
831  for (int i = 0; i < p; i++) // (2,3)
832  {
833  Nodes.IntPoint(o).Set3(0., eop[pm1-i], eop[i]);
834  dof2tk[o++] = 5;
835  }
836 
837  // faces
838  for (int j = 0; j <= pm2; j++) // (1,2,3)
839  for (int i = 0; i + j <= pm2; i++)
840  {
841  double w = fop[i] + fop[j] + fop[pm2-i-j];
842  Nodes.IntPoint(o).Set3(fop[pm2-i-j]/w, fop[i]/w, fop[j]/w);
843  dof2tk[o++] = 3;
844  Nodes.IntPoint(o).Set3(fop[pm2-i-j]/w, fop[i]/w, fop[j]/w);
845  dof2tk[o++] = 4;
846  }
847  for (int j = 0; j <= pm2; j++) // (0,3,2)
848  for (int i = 0; i + j <= pm2; i++)
849  {
850  double w = fop[i] + fop[j] + fop[pm2-i-j];
851  Nodes.IntPoint(o).Set3(0., fop[j]/w, fop[i]/w);
852  dof2tk[o++] = 2;
853  Nodes.IntPoint(o).Set3(0., fop[j]/w, fop[i]/w);
854  dof2tk[o++] = 1;
855  }
856  for (int j = 0; j <= pm2; j++) // (0,1,3)
857  for (int i = 0; i + j <= pm2; i++)
858  {
859  double w = fop[i] + fop[j] + fop[pm2-i-j];
860  Nodes.IntPoint(o).Set3(fop[i]/w, 0., fop[j]/w);
861  dof2tk[o++] = 0;
862  Nodes.IntPoint(o).Set3(fop[i]/w, 0., fop[j]/w);
863  dof2tk[o++] = 2;
864  }
865  for (int j = 0; j <= pm2; j++) // (0,2,1)
866  for (int i = 0; i + j <= pm2; i++)
867  {
868  double w = fop[i] + fop[j] + fop[pm2-i-j];
869  Nodes.IntPoint(o).Set3(fop[j]/w, fop[i]/w, 0.);
870  dof2tk[o++] = 1;
871  Nodes.IntPoint(o).Set3(fop[j]/w, fop[i]/w, 0.);
872  dof2tk[o++] = 0;
873  }
874 
875  // interior
876  for (int k = 0; k <= pm3; k++)
877  for (int j = 0; j + k <= pm3; j++)
878  for (int i = 0; i + j + k <= pm3; i++)
879  {
880  double w = iop[i] + iop[j] + iop[k] + iop[pm3-i-j-k];
881  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
882  dof2tk[o++] = 0;
883  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
884  dof2tk[o++] = 1;
885  Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
886  dof2tk[o++] = 2;
887  }
888 
889  DenseMatrix T(dof);
890  for (int m = 0; m < dof; m++)
891  {
892  const IntegrationPoint &ip = Nodes.IntPoint(m);
893  const double *tm = tk + 3*dof2tk[m];
894  o = 0;
895 
896  poly1d.CalcBasis(pm1, ip.x, shape_x);
897  poly1d.CalcBasis(pm1, ip.y, shape_y);
898  poly1d.CalcBasis(pm1, ip.z, shape_z);
899  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y - ip.z, shape_l);
900 
901  for (int k = 0; k <= pm1; k++)
902  for (int j = 0; j + k <= pm1; j++)
903  for (int i = 0; i + j + k <= pm1; i++)
904  {
905  double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(pm1-i-j-k);
906  T(o++, m) = s * tm[0];
907  T(o++, m) = s * tm[1];
908  T(o++, m) = s * tm[2];
909  }
910  for (int k = 0; k <= pm1; k++)
911  for (int j = 0; j + k <= pm1; j++)
912  {
913  double s = shape_x(pm1-j-k)*shape_y(j)*shape_z(k);
914  T(o++, m) = s*((ip.y - c)*tm[0] - (ip.x - c)*tm[1]);
915  T(o++, m) = s*((ip.z - c)*tm[0] - (ip.x - c)*tm[2]);
916  }
917  for (int k = 0; k <= pm1; k++)
918  {
919  T(o++, m) =
920  shape_y(pm1-k)*shape_z(k)*((ip.z - c)*tm[1] - (ip.y - c)*tm[2]);
921  }
922  }
923 
924  Ti.Factor(T);
925  // mfem::out << "ND_TetrahedronElement(" << p << ") : "; Ti.TestInversion();
926 }
927 
929  DenseMatrix &shape) const
930 {
931  const int pm1 = order - 1;
932 
933 #ifdef MFEM_THREAD_SAFE
934  const int p = order;
935  Vector shape_x(p), shape_y(p), shape_z(p), shape_l(p);
936  DenseMatrix u(dof, dim);
937 #endif
938 
939  poly1d.CalcBasis(pm1, ip.x, shape_x);
940  poly1d.CalcBasis(pm1, ip.y, shape_y);
941  poly1d.CalcBasis(pm1, ip.z, shape_z);
942  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y - ip.z, shape_l);
943 
944  int n = 0;
945  for (int k = 0; k <= pm1; k++)
946  for (int j = 0; j + k <= pm1; j++)
947  for (int i = 0; i + j + k <= pm1; i++)
948  {
949  double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(pm1-i-j-k);
950  u(n,0) = s; u(n,1) = 0.; u(n,2) = 0.; n++;
951  u(n,0) = 0.; u(n,1) = s; u(n,2) = 0.; n++;
952  u(n,0) = 0.; u(n,1) = 0.; u(n,2) = s; n++;
953  }
954  for (int k = 0; k <= pm1; k++)
955  for (int j = 0; j + k <= pm1; j++)
956  {
957  double s = shape_x(pm1-j-k)*shape_y(j)*shape_z(k);
958  u(n,0) = s*(ip.y - c); u(n,1) = -s*(ip.x - c); u(n,2) = 0.; n++;
959  u(n,0) = s*(ip.z - c); u(n,1) = 0.; u(n,2) = -s*(ip.x - c); n++;
960  }
961  for (int k = 0; k <= pm1; k++)
962  {
963  double s = shape_y(pm1-k)*shape_z(k);
964  u(n,0) = 0.; u(n,1) = s*(ip.z - c); u(n,2) = -s*(ip.y - c); n++;
965  }
966 
967  Ti.Mult(u, shape);
968 }
969 
971  DenseMatrix &curl_shape) const
972 {
973  const int pm1 = order - 1;
974 
975 #ifdef MFEM_THREAD_SAFE
976  const int p = order;
977  Vector shape_x(p), shape_y(p), shape_z(p), shape_l(p);
978  Vector dshape_x(p), dshape_y(p), dshape_z(p), dshape_l(p);
979  DenseMatrix u(dof, dim);
980 #endif
981 
982  poly1d.CalcBasis(pm1, ip.x, shape_x, dshape_x);
983  poly1d.CalcBasis(pm1, ip.y, shape_y, dshape_y);
984  poly1d.CalcBasis(pm1, ip.z, shape_z, dshape_z);
985  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y - ip.z, shape_l, dshape_l);
986 
987  int n = 0;
988  for (int k = 0; k <= pm1; k++)
989  for (int j = 0; j + k <= pm1; j++)
990  for (int i = 0; i + j + k <= pm1; i++)
991  {
992  int l = pm1-i-j-k;
993  const double dx = (dshape_x(i)*shape_l(l) -
994  shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
995  const double dy = (dshape_y(j)*shape_l(l) -
996  shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
997  const double dz = (dshape_z(k)*shape_l(l) -
998  shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
999 
1000  u(n,0) = 0.; u(n,1) = dz; u(n,2) = -dy; n++;
1001  u(n,0) = -dz; u(n,1) = 0.; u(n,2) = dx; n++;
1002  u(n,0) = dy; u(n,1) = -dx; u(n,2) = 0.; n++;
1003  }
1004  for (int k = 0; k <= pm1; k++)
1005  for (int j = 0; j + k <= pm1; j++)
1006  {
1007  int i = pm1 - j - k;
1008  // s = shape_x(i)*shape_y(j)*shape_z(k);
1009  // curl of s*(ip.y - c, -(ip.x - c), 0):
1010  u(n,0) = shape_x(i)*(ip.x - c)*shape_y(j)*dshape_z(k);
1011  u(n,1) = shape_x(i)*shape_y(j)*(ip.y - c)*dshape_z(k);
1012  u(n,2) =
1013  -((dshape_x(i)*(ip.x - c) + shape_x(i))*shape_y(j)*shape_z(k) +
1014  (dshape_y(j)*(ip.y - c) + shape_y(j))*shape_x(i)*shape_z(k));
1015  n++;
1016  // curl of s*(ip.z - c, 0, -(ip.x - c)):
1017  u(n,0) = -shape_x(i)*(ip.x - c)*dshape_y(j)*shape_z(k);
1018  u(n,1) = (shape_x(i)*shape_y(j)*(dshape_z(k)*(ip.z - c) + shape_z(k)) +
1019  (dshape_x(i)*(ip.x - c) + shape_x(i))*shape_y(j)*shape_z(k));
1020  u(n,2) = -shape_x(i)*dshape_y(j)*shape_z(k)*(ip.z - c);
1021  n++;
1022  }
1023  for (int k = 0; k <= pm1; k++)
1024  {
1025  int j = pm1 - k;
1026  // curl of shape_y(j)*shape_z(k)*(0, ip.z - c, -(ip.y - c)):
1027  u(n,0) = -((dshape_y(j)*(ip.y - c) + shape_y(j))*shape_z(k) +
1028  shape_y(j)*(dshape_z(k)*(ip.z - c) + shape_z(k)));
1029  u(n,1) = 0.;
1030  u(n,2) = 0.; n++;
1031  }
1032 
1033  Ti.Mult(u, curl_shape);
1034 }
1035 
1036 
1037 const double ND_TriangleElement::tk[8] =
1038 { 1.,0., -1.,1., 0.,-1., 0.,1. };
1039 
1040 const double ND_TriangleElement::c = 1./3.;
1041 
1043  : VectorFiniteElement(2, Geometry::TRIANGLE, p*(p + 2), p,
1044  H_CURL, FunctionSpace::Pk),
1045  dof2tk(dof)
1046 {
1047  const double *eop = poly1d.OpenPoints(p - 1);
1048  const double *iop = (p > 1) ? poly1d.OpenPoints(p - 2) : NULL;
1049 
1050  const int pm1 = p - 1, pm2 = p - 2;
1051 
1052 #ifndef MFEM_THREAD_SAFE
1053  shape_x.SetSize(p);
1054  shape_y.SetSize(p);
1055  shape_l.SetSize(p);
1056  dshape_x.SetSize(p);
1057  dshape_y.SetSize(p);
1058  dshape_l.SetSize(p);
1059  u.SetSize(dof, dim);
1060  curlu.SetSize(dof);
1061 #else
1062  Vector shape_x(p), shape_y(p), shape_l(p);
1063 #endif
1064 
1065  int n = 0;
1066  // edges
1067  for (int i = 0; i < p; i++) // (0,1)
1068  {
1069  Nodes.IntPoint(n).Set2(eop[i], 0.);
1070  dof2tk[n++] = 0;
1071  }
1072  for (int i = 0; i < p; i++) // (1,2)
1073  {
1074  Nodes.IntPoint(n).Set2(eop[pm1-i], eop[i]);
1075  dof2tk[n++] = 1;
1076  }
1077  for (int i = 0; i < p; i++) // (2,0)
1078  {
1079  Nodes.IntPoint(n).Set2(0., eop[pm1-i]);
1080  dof2tk[n++] = 2;
1081  }
1082 
1083  // interior
1084  for (int j = 0; j <= pm2; j++)
1085  for (int i = 0; i + j <= pm2; i++)
1086  {
1087  double w = iop[i] + iop[j] + iop[pm2-i-j];
1088  Nodes.IntPoint(n).Set2(iop[i]/w, iop[j]/w);
1089  dof2tk[n++] = 0;
1090  Nodes.IntPoint(n).Set2(iop[i]/w, iop[j]/w);
1091  dof2tk[n++] = 3;
1092  }
1093 
1094  DenseMatrix T(dof);
1095  for (int m = 0; m < dof; m++)
1096  {
1097  const IntegrationPoint &ip = Nodes.IntPoint(m);
1098  const double *tm = tk + 2*dof2tk[m];
1099  n = 0;
1100 
1101  poly1d.CalcBasis(pm1, ip.x, shape_x);
1102  poly1d.CalcBasis(pm1, ip.y, shape_y);
1103  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y, shape_l);
1104 
1105  for (int j = 0; j <= pm1; j++)
1106  for (int i = 0; i + j <= pm1; i++)
1107  {
1108  double s = shape_x(i)*shape_y(j)*shape_l(pm1-i-j);
1109  T(n++, m) = s * tm[0];
1110  T(n++, m) = s * tm[1];
1111  }
1112  for (int j = 0; j <= pm1; j++)
1113  {
1114  T(n++, m) =
1115  shape_x(pm1-j)*shape_y(j)*((ip.y - c)*tm[0] - (ip.x - c)*tm[1]);
1116  }
1117  }
1118 
1119  Ti.Factor(T);
1120  // mfem::out << "ND_TriangleElement(" << p << ") : "; Ti.TestInversion();
1121 }
1122 
1124  DenseMatrix &shape) const
1125 {
1126  const int pm1 = order - 1;
1127 
1128 #ifdef MFEM_THREAD_SAFE
1129  const int p = order;
1130  Vector shape_x(p), shape_y(p), shape_l(p);
1131  DenseMatrix u(dof, dim);
1132 #endif
1133 
1134  poly1d.CalcBasis(pm1, ip.x, shape_x);
1135  poly1d.CalcBasis(pm1, ip.y, shape_y);
1136  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y, shape_l);
1137 
1138  int n = 0;
1139  for (int j = 0; j <= pm1; j++)
1140  for (int i = 0; i + j <= pm1; i++)
1141  {
1142  double s = shape_x(i)*shape_y(j)*shape_l(pm1-i-j);
1143  u(n,0) = s; u(n,1) = 0; n++;
1144  u(n,0) = 0; u(n,1) = s; n++;
1145  }
1146  for (int j = 0; j <= pm1; j++)
1147  {
1148  double s = shape_x(pm1-j)*shape_y(j);
1149  u(n,0) = s*(ip.y - c);
1150  u(n,1) = -s*(ip.x - c);
1151  n++;
1152  }
1153 
1154  Ti.Mult(u, shape);
1155 }
1156 
1158  DenseMatrix &curl_shape) const
1159 {
1160  const int pm1 = order - 1;
1161 
1162 #ifdef MFEM_THREAD_SAFE
1163  const int p = order;
1164  Vector shape_x(p), shape_y(p), shape_l(p);
1165  Vector dshape_x(p), dshape_y(p), dshape_l(p);
1166  Vector curlu(dof);
1167 #endif
1168 
1169  poly1d.CalcBasis(pm1, ip.x, shape_x, dshape_x);
1170  poly1d.CalcBasis(pm1, ip.y, shape_y, dshape_y);
1171  poly1d.CalcBasis(pm1, 1. - ip.x - ip.y, shape_l, dshape_l);
1172 
1173  int n = 0;
1174  for (int j = 0; j <= pm1; j++)
1175  for (int i = 0; i + j <= pm1; i++)
1176  {
1177  int l = pm1-i-j;
1178  const double dx = (dshape_x(i)*shape_l(l) -
1179  shape_x(i)*dshape_l(l)) * shape_y(j);
1180  const double dy = (dshape_y(j)*shape_l(l) -
1181  shape_y(j)*dshape_l(l)) * shape_x(i);
1182 
1183  curlu(n++) = -dy;
1184  curlu(n++) = dx;
1185  }
1186 
1187  for (int j = 0; j <= pm1; j++)
1188  {
1189  int i = pm1 - j;
1190  // curl of shape_x(i)*shape_y(j) * (ip.y - c, -(ip.x - c), 0):
1191  curlu(n++) = -((dshape_x(i)*(ip.x - c) + shape_x(i)) * shape_y(j) +
1192  (dshape_y(j)*(ip.y - c) + shape_y(j)) * shape_x(i));
1193  }
1194 
1195  Vector curl2d(curl_shape.Data(),dof);
1196  Ti.Mult(curlu, curl2d);
1197 }
1198 
1199 
1200 const double ND_SegmentElement::tk[1] = { 1. };
1201 
1202 ND_SegmentElement::ND_SegmentElement(const int p, const int ob_type)
1203  : VectorTensorFiniteElement(1, p, p - 1, ob_type, H_CURL,
1204  DofMapType::L2_DOF_MAP),
1205  dof2tk(dof)
1206 {
1207  if (obasis1d.IsIntegratedType()) { is_nodal = false; }
1208 
1209  const double *op = poly1d.OpenPoints(p - 1, ob_type);
1210 
1211  // set dof2tk and Nodes
1212  for (int i = 0; i < p; i++)
1213  {
1214  dof2tk[i] = 0;
1215  Nodes.IntPoint(i).x = op[i];
1216  }
1217 }
1218 
1220  DenseMatrix &shape) const
1221 {
1222  Vector vshape(shape.Data(), dof);
1223 
1224  obasis1d.Eval(ip.x, vshape);
1225 }
1226 
1227 const double ND_WedgeElement::tk[15] =
1228 { 1.,0.,0., -1.,1.,0., 0.,-1.,0., 0.,0.,1., 0.,1.,0. };
1229 
1231  const int cb_type,
1232  const int ob_type)
1233  : VectorFiniteElement(3, Geometry::PRISM,
1234  3 * p * ((p + 1) * (p + 2))/2, p,
1235  H_CURL, FunctionSpace::Qk),
1236  dof2tk(dof),
1237  t_dof(dof),
1238  s_dof(dof),
1239  H1TriangleFE(p, cb_type),
1240  NDTriangleFE(p),
1241  H1SegmentFE(p, cb_type),
1242  NDSegmentFE(p, ob_type)
1243 {
1244  MFEM_ASSERT(H1TriangleFE.GetDof() * NDSegmentFE.GetDof() +
1245  NDTriangleFE.GetDof() * H1SegmentFE.GetDof() == dof,
1246  "Mismatch in number of degrees of freedom "
1247  "when building ND_WedgeElement!");
1248 
1249 #ifndef MFEM_THREAD_SAFE
1250  t1_shape.SetSize(H1TriangleFE.GetDof());
1251  s1_shape.SetSize(H1SegmentFE.GetDof());
1252  tn_shape.SetSize(NDTriangleFE.GetDof(), 2);
1253  sn_shape.SetSize(NDSegmentFE.GetDof(), 1);
1254  t1_dshape.SetSize(H1TriangleFE.GetDof(), 2);
1255  s1_dshape.SetSize(H1SegmentFE.GetDof(), 1);
1256  tn_dshape.SetSize(NDTriangleFE.GetDof(), 1);
1257 #endif
1258 
1259  const int pm1 = p - 1, pm2 = p - 2;
1260 
1261  const IntegrationRule &t1_n = H1TriangleFE.GetNodes();
1262  const IntegrationRule &tn_n = NDTriangleFE.GetNodes();
1263  const IntegrationRule &s1_n = H1SegmentFE.GetNodes();
1264  const IntegrationRule &sn_n = NDSegmentFE.GetNodes();
1265 
1266  // edges
1267  int o = 0;
1268  for (int i = 0; i < p; i++) // (0,1)
1269  {
1270  t_dof[o] = i; s_dof[o] = 0; dof2tk[o] = 0;
1271  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1272  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1273  o++;
1274  }
1275  for (int i = 0; i < p; i++) // (1,2)
1276  {
1277  t_dof[o] = p + i; s_dof[o] = 0; dof2tk[o] = 1;
1278  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1279  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1280  o++;
1281  }
1282  for (int i = 0; i < p; i++) // (2,0)
1283  {
1284  t_dof[o] = 2 * p + i; s_dof[o] = 0; dof2tk[o] = 2;
1285  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1286  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1287  o++;
1288  }
1289  for (int i = 0; i < p; i++) // (3,4)
1290  {
1291  t_dof[o] = i; s_dof[o] = 1; dof2tk[o] = 0;
1292  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1293  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1294  o++;
1295  }
1296  for (int i = 0; i < p; i++) // (4,5)
1297  {
1298  t_dof[o] = p + i; s_dof[o] = 1; dof2tk[o] = 1;
1299  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1300  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1301  o++;
1302  }
1303  for (int i = 0; i < p; i++) // (5,3)
1304  {
1305  t_dof[o] = 2 * p + i; s_dof[o] = 1; dof2tk[o] = 2;
1306  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1307  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1308  o++;
1309  }
1310  for (int i = 0; i < p; i++) // (0,3)
1311  {
1312  t_dof[o] = 0; s_dof[o] = i; dof2tk[o] = 3;
1313  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1314  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1315  o++;
1316  }
1317  for (int i = 0; i < p; i++) // (1,4)
1318  {
1319  t_dof[o] = 1; s_dof[o] = i; dof2tk[o] = 3;
1320  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1321  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1322  o++;
1323  }
1324  for (int i = 0; i < p; i++) // (2,5)
1325  {
1326  t_dof[o] = 2; s_dof[o] = i; dof2tk[o] = 3;
1327  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1328  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1329  o++;
1330  }
1331 
1332  // faces
1333  // (0,2,1) -- bottom
1334  int l = 0;
1335  for (int j = 0; j <= pm2; j++)
1336  for (int i = 0; i + j <= pm2; i++)
1337  {
1338  l = j + ( 2 * p - 1 - i) * i / 2;
1339  t_dof[o] = 3 * p + 2*l+1; s_dof[o] = 0; dof2tk[o] = 4;
1340  const IntegrationPoint & t_ip0 = tn_n.IntPoint(t_dof[o]);
1341  Nodes.IntPoint(o).Set3(t_ip0.x, t_ip0.y, s1_n.IntPoint(s_dof[o]).x);
1342  o++;
1343  t_dof[o] = 3 * p + 2*l; s_dof[o] = 0; dof2tk[o] = 0;
1344  const IntegrationPoint & t_ip1 = tn_n.IntPoint(t_dof[o]);
1345  Nodes.IntPoint(o).Set3(t_ip1.x, t_ip1.y, s1_n.IntPoint(s_dof[o]).x);
1346  o++;
1347  }
1348  // (3,4,5) -- top
1349  int m = 0;
1350  for (int j = 0; j <= pm2; j++)
1351  for (int i = 0; i + j <= pm2; i++)
1352  {
1353  t_dof[o] = 3 * p + m; s_dof[o] = 1; dof2tk[o] = 0; m++;
1354  const IntegrationPoint & t_ip0 = tn_n.IntPoint(t_dof[o]);
1355  Nodes.IntPoint(o).Set3(t_ip0.x, t_ip0.y, s1_n.IntPoint(s_dof[o]).x);
1356  o++;
1357  t_dof[o] = 3 * p + m; s_dof[o] = 1; dof2tk[o] = 4; m++;
1358  const IntegrationPoint & t_ip1 = tn_n.IntPoint(t_dof[o]);
1359  Nodes.IntPoint(o).Set3(t_ip1.x, t_ip1.y, s1_n.IntPoint(s_dof[o]).x);
1360  o++;
1361  }
1362  // (0, 1, 4, 3) -- xz plane
1363  for (int j = 2; j <= p; j++)
1364  for (int i = 0; i < p; i++)
1365  {
1366  t_dof[o] = i; s_dof[o] = j; dof2tk[o] = 0;
1367  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1368  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1369  o++;
1370  }
1371  for (int j = 0; j < p; j++)
1372  for (int i = 0; i < pm1; i++)
1373  {
1374  t_dof[o] = 3 + i; s_dof[o] = j; dof2tk[o] = 3;
1375  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1376  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1377  o++;
1378  }
1379  // (1, 2, 5, 4) -- (y-x)z plane
1380  for (int j = 2; j <= p; j++)
1381  for (int i = 0; i < p; i++)
1382  {
1383  t_dof[o] = p + i; s_dof[o] = j; dof2tk[o] = 1;
1384  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1385  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1386  o++;
1387  }
1388  for (int j = 0; j < p; j++)
1389  for (int i = 0; i < pm1; i++)
1390  {
1391  t_dof[o] = p + 2 + i; s_dof[o] = j; dof2tk[o] = 3;
1392  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1393  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1394  o++;
1395  }
1396  // (2, 0, 3, 5) -- yz plane
1397  for (int j = 2; j <= p; j++)
1398  for (int i = 0; i < p; i++)
1399  {
1400  t_dof[o] = 2 * p + i; s_dof[o] = j; dof2tk[o] = 2;
1401  const IntegrationPoint & t_ip = tn_n.IntPoint(t_dof[o]);
1402  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, s1_n.IntPoint(s_dof[o]).x);
1403  o++;
1404  }
1405  for (int j = 0; j < p; j++)
1406  for (int i = 0; i < pm1; i++)
1407  {
1408  t_dof[o] = 2 * p + 1 + i; s_dof[o] = j; dof2tk[o] = 3;
1409  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1410  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1411  o++;
1412  }
1413 
1414  // interior
1415  for (int k = 2; k <= p; k++)
1416  {
1417  l = 0;
1418  for (int j = 0; j <= pm2; j++)
1419  for (int i = 0; i + j <= pm2; i++)
1420  {
1421  t_dof[o] = 3 * p + l; s_dof[o] = k; dof2tk[o] = 0; l++;
1422  const IntegrationPoint & t_ip0 = tn_n.IntPoint(t_dof[o]);
1423  Nodes.IntPoint(o).Set3(t_ip0.x, t_ip0.y, s1_n.IntPoint(s_dof[o]).x);
1424  o++;
1425  t_dof[o] = 3 * p + l; s_dof[o] = k; dof2tk[o] = 4; l++;
1426  const IntegrationPoint & t_ip1 = tn_n.IntPoint(t_dof[o]);
1427  Nodes.IntPoint(o).Set3(t_ip1.x, t_ip1.y, s1_n.IntPoint(s_dof[o]).x);
1428  o++;
1429  }
1430  }
1431  for (int k = 0; k < p; k++)
1432  {
1433  l = 0;
1434  for (int j = 0; j < pm2; j++)
1435  for (int i = 0; i + j < pm2; i++)
1436  {
1437  t_dof[o] = 3 * p + l; s_dof[o] = k; dof2tk[o] = 3; l++;
1438  const IntegrationPoint & t_ip = t1_n.IntPoint(t_dof[o]);
1439  Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sn_n.IntPoint(s_dof[o]).x);
1440  o++;
1441  }
1442  }
1443 }
1444 
1446  DenseMatrix &shape) const
1447 {
1448 #ifdef MFEM_THREAD_SAFE
1449  Vector t1_shape(H1TriangleFE.GetDof());
1450  Vector s1_shape(H1SegmentFE.GetDof());
1451  DenseMatrix tn_shape(NDTriangleFE.GetDof(), 2);
1452  DenseMatrix sn_shape(NDSegmentFE.GetDof(), 1);
1453 #endif
1454 
1455  IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1456 
1457  H1TriangleFE.CalcShape(ip, t1_shape);
1458  NDTriangleFE.CalcVShape(ip, tn_shape);
1459  H1SegmentFE.CalcShape(ipz, s1_shape);
1460  NDSegmentFE.CalcVShape(ipz, sn_shape);
1461 
1462  for (int i=0; i<dof; i++)
1463  {
1464  if ( dof2tk[i] != 3 )
1465  {
1466  shape(i, 0) = tn_shape(t_dof[i], 0) * s1_shape[s_dof[i]];
1467  shape(i, 1) = tn_shape(t_dof[i], 1) * s1_shape[s_dof[i]];
1468  shape(i, 2) = 0.0;
1469  }
1470  else
1471  {
1472  shape(i, 0) = 0.0;
1473  shape(i, 1) = 0.0;
1474  shape(i, 2) = t1_shape[t_dof[i]] * sn_shape(s_dof[i], 0);
1475  }
1476  }
1477 }
1478 
1480  DenseMatrix &curl_shape) const
1481 {
1482 #ifdef MFEM_THREAD_SAFE
1483  Vector s1_shape(H1SegmentFE.GetDof());
1484  DenseMatrix t1_dshape(H1TriangleFE.GetDof(), 2);
1485  DenseMatrix s1_dshape(H1SegmentFE.GetDof(), 1);
1486  DenseMatrix tn_shape(NDTriangleFE.GetDof(), 2);
1487  DenseMatrix sn_shape(NDSegmentFE.GetDof(), 1);
1488  DenseMatrix tn_dshape(NDTriangleFE.GetDof(), 1);
1489 #endif
1490 
1491  IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1492 
1493  H1TriangleFE.CalcDShape(ip, t1_dshape);
1494  H1SegmentFE.CalcShape(ipz, s1_shape);
1495  H1SegmentFE.CalcDShape(ipz, s1_dshape);
1496  NDTriangleFE.CalcVShape(ip, tn_shape);
1497  NDTriangleFE.CalcCurlShape(ip, tn_dshape);
1498  NDSegmentFE.CalcVShape(ipz, sn_shape);
1499 
1500  for (int i=0; i<dof; i++)
1501  {
1502  if ( dof2tk[i] != 3 )
1503  {
1504  curl_shape(i, 0) = -tn_shape(t_dof[i], 1) * s1_dshape(s_dof[i], 0);
1505  curl_shape(i, 1) = tn_shape(t_dof[i], 0) * s1_dshape(s_dof[i], 0);
1506  curl_shape(i, 2) = tn_dshape(t_dof[i], 0) * s1_shape[s_dof[i]];
1507  }
1508  else
1509  {
1510  curl_shape(i, 0) = t1_dshape(t_dof[i], 1) * sn_shape(s_dof[i], 0);
1511  curl_shape(i, 1) = -t1_dshape(t_dof[i], 0) * sn_shape(s_dof[i], 0);
1512  curl_shape(i, 2) = 0.0;
1513  }
1514  }
1515 }
1516 
1518  : VectorFiniteElement(1, Geometry::POINT, 2, p,
1519  H_CURL, FunctionSpace::Pk)
1520 {
1521  // VectorFiniteElement::SetDerivMembers doesn't support 0D H_CURL elements
1522  // so we mimic a 1D element and then correct the dimension here.
1523  dim = 0;
1524  vdim = 2;
1525  cdim = 0;
1526 }
1527 
1529  DenseMatrix &shape) const
1530 {
1531  shape(0,0) = 1.0;
1532  shape(0,1) = 0.0;
1533 
1534  shape(1,0) = 0.0;
1535  shape(1,1) = 1.0;
1536 }
1537 
1539  DenseMatrix &shape) const
1540 {
1541  CalcVShape(Trans.GetIntPoint(), shape);
1542 }
1543 
1544 const double ND_R1D_SegmentElement::tk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1545 
1547  const int cb_type,
1548  const int ob_type)
1549  : VectorFiniteElement(1, Geometry::SEGMENT, 3 * p + 2, p,
1550  H_CURL, FunctionSpace::Pk),
1551  dof2tk(dof),
1552  cbasis1d(poly1d.GetBasis(p, VerifyClosed(cb_type))),
1553  obasis1d(poly1d.GetBasis(p - 1, VerifyOpen(ob_type)))
1554 {
1555  // Override default types for VectorFiniteElements
1556  deriv_type = CURL;
1559 
1560  // Override default dimensions for VectorFiniteElements
1561  vdim = 3;
1562  cdim = 3;
1563 
1564  const double *cp = poly1d.ClosedPoints(p, cb_type);
1565  const double *op = poly1d.OpenPoints(p - 1, ob_type);
1566 
1567 #ifndef MFEM_THREAD_SAFE
1568  shape_cx.SetSize(p + 1);
1569  shape_ox.SetSize(p);
1570  dshape_cx.SetSize(p + 1);
1571 #endif
1572 
1573  dof_map.SetSize(dof);
1574 
1575  int o = 0;
1576  // nodes
1577  // (0)
1578  Nodes.IntPoint(o).x = cp[0]; // y-directed
1579  dof_map[p] = o; dof2tk[o++] = 1;
1580  Nodes.IntPoint(o).x = cp[0]; // z-directed
1581  dof_map[2*p+1] = o; dof2tk[o++] = 2;
1582 
1583  // (1)
1584  Nodes.IntPoint(o).x = cp[p]; // y-directed
1585  dof_map[2*p] = o; dof2tk[o++] = 1;
1586  Nodes.IntPoint(o).x = cp[p]; // z-directed
1587  dof_map[3*p+1] = o; dof2tk[o++] = 2;
1588 
1589  // interior
1590  // x-components
1591  for (int i = 0; i < p; i++)
1592  {
1593  Nodes.IntPoint(o).x = op[i];
1594  dof_map[i] = o; dof2tk[o++] = 0;
1595  }
1596  // y-components
1597  for (int i = 1; i < p; i++)
1598  {
1599  Nodes.IntPoint(o).x = cp[i];
1600  dof_map[p+i] = o; dof2tk[o++] = 1;
1601  }
1602  // z-components
1603  for (int i = 1; i < p; i++)
1604  {
1605  Nodes.IntPoint(o).x = cp[i];
1606  dof_map[2*p+1+i] = o; dof2tk[o++] = 2;
1607  }
1608 }
1609 
1611  DenseMatrix &shape) const
1612 {
1613  const int p = order;
1614 
1615 #ifdef MFEM_THREAD_SAFE
1616  Vector shape_cx(p + 1), shape_ox(p);
1617 #endif
1618 
1619  cbasis1d.Eval(ip.x, shape_cx);
1620  obasis1d.Eval(ip.x, shape_ox);
1621 
1622  int o = 0;
1623  // x-components
1624  for (int i = 0; i < p; i++)
1625  {
1626  int idx = dof_map[o++];
1627  shape(idx,0) = shape_ox(i);
1628  shape(idx,1) = 0.;
1629  shape(idx,2) = 0.;
1630  }
1631  // y-components
1632  for (int i = 0; i <= p; i++)
1633  {
1634  int idx = dof_map[o++];
1635  shape(idx,0) = 0.;
1636  shape(idx,1) = shape_cx(i);
1637  shape(idx,2) = 0.;
1638  }
1639  // z-components
1640  for (int i = 0; i <= p; i++)
1641  {
1642  int idx = dof_map[o++];
1643  shape(idx,0) = 0.;
1644  shape(idx,1) = 0.;
1645  shape(idx,2) = shape_cx(i);
1646  }
1647 }
1648 
1650  DenseMatrix &shape) const
1651 {
1652  CalcVShape(Trans.GetIntPoint(), shape);
1653  const DenseMatrix & JI = Trans.InverseJacobian();
1654  MFEM_ASSERT(JI.Width() == 1 && JI.Height() == 1,
1655  "ND_R1D_SegmentElement cannot be embedded in "
1656  "2 or 3 dimensional spaces");
1657  for (int i=0; i<dof; i++)
1658  {
1659  shape(i, 0) *= JI(0,0);
1660  }
1661 }
1662 
1664  DenseMatrix &curl_shape) const
1665 {
1666  const int p = order;
1667 
1668 #ifdef MFEM_THREAD_SAFE
1669  Vector shape_cx(p + 1), shape_ox(p);
1670  Vector dshape_cx(p + 1);
1671 #endif
1672 
1673  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
1674  obasis1d.Eval(ip.x, shape_ox);
1675 
1676  int o = 0;
1677  // x-components
1678  for (int i = 0; i < p; i++)
1679  {
1680  int idx = dof_map[o++];
1681  curl_shape(idx,0) = 0.;
1682  curl_shape(idx,1) = 0.;
1683  curl_shape(idx,2) = 0.;
1684  }
1685  // y-components
1686  for (int i = 0; i <= p; i++)
1687  {
1688  int idx = dof_map[o++];
1689  curl_shape(idx,0) = 0.;
1690  curl_shape(idx,1) = 0.;
1691  curl_shape(idx,2) = dshape_cx(i);
1692  }
1693  // z-components
1694  for (int i = 0; i <= p; i++)
1695  {
1696  int idx = dof_map[o++];
1697  curl_shape(idx,0) = 0.;
1698  curl_shape(idx,1) = -dshape_cx(i);
1699  curl_shape(idx,2) = 0.;
1700  }
1701 }
1702 
1704  DenseMatrix &curl_shape) const
1705 {
1706  CalcCurlShape(Trans.GetIntPoint(), curl_shape);
1707  const DenseMatrix & J = Trans.Jacobian();
1708  MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1709  "ND_R1D_SegmentElement cannot be embedded in "
1710  "2 or 3 dimensional spaces");
1711  curl_shape *= (1.0 / J.Weight());
1712 }
1713 
1716  Vector &dofs) const
1717 {
1718  double data[3];
1719  Vector vk(data, 3);
1720 
1721  for (int k = 0; k < dof; k++)
1722  {
1723  Trans.SetIntPoint(&Nodes.IntPoint(k));
1724 
1725  vc.Eval(vk, Trans, Nodes.IntPoint(k));
1726  // dof_k = vk^t J tk
1727  Vector t(const_cast<double*>(&tk[dof2tk[k] * 3]), 3);
1728  dofs(k) = Trans.Jacobian()(0,0) * t(0) * vk(0) +
1729  t(1) * vk(1) + t(2) * vk(2);
1730  }
1731 
1732 }
1733 
1736  DenseMatrix &I) const
1737 {
1738  if (fe.GetRangeType() == SCALAR)
1739  {
1740  double vk[Geometry::MaxDim];
1741  Vector shape(fe.GetDof());
1742 
1743  double * tk_ptr = const_cast<double*>(tk);
1744 
1745  I.SetSize(dof, vdim*fe.GetDof());
1746  for (int k = 0; k < dof; k++)
1747  {
1748  const IntegrationPoint &ip = Nodes.IntPoint(k);
1749 
1750  Vector t1(&tk_ptr[dof2tk[k] * 3], 1);
1751  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
1752 
1753  fe.CalcShape(ip, shape);
1754  Trans.SetIntPoint(&ip);
1755  // Transform ND edge tengents from reference to physical space
1756  // vk = J tk
1757  Trans.Jacobian().Mult(t1, vk);
1758  vk[1] = t3[1];
1759  vk[2] = t3[2];
1760  if (fe.GetMapType() == INTEGRAL)
1761  {
1762  double w = 1.0/Trans.Weight();
1763  for (int d = 0; d < vdim; d++)
1764  {
1765  vk[d] *= w;
1766  }
1767  }
1768 
1769  for (int j = 0; j < shape.Size(); j++)
1770  {
1771  double s = shape(j);
1772  if (fabs(s) < 1e-12)
1773  {
1774  s = 0.0;
1775  }
1776  // Project scalar basis function multiplied by each coordinate
1777  // direction onto the transformed edge tangents
1778  for (int d = 0; d < vdim; d++)
1779  {
1780  I(k, j + d*shape.Size()) = s*vk[d];
1781  }
1782  }
1783  }
1784  }
1785  else
1786  {
1787  double vk[Geometry::MaxDim];
1788  DenseMatrix vshape(fe.GetDof(), fe.GetVDim());
1789 
1790  double * tk_ptr = const_cast<double*>(tk);
1791 
1792  I.SetSize(dof, fe.GetDof());
1793  for (int k = 0; k < dof; k++)
1794  {
1795  const IntegrationPoint &ip = Nodes.IntPoint(k);
1796 
1797  Vector t1(&tk_ptr[dof2tk[k] * 3], 1);
1798  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
1799 
1800  Trans.SetIntPoint(&ip);
1801  // Transform ND edge tangents from reference to physical space
1802  // vk = J tk
1803  Trans.Jacobian().Mult(t1, vk);
1804  // Compute fe basis functions in physical space
1805  fe.CalcVShape(Trans, vshape);
1806  // Project fe basis functions onto transformed edge tangents
1807  for (int j=0; j<vshape.Height(); j++)
1808  {
1809  I(k, j) = 0.0;
1810  I(k, j) += vshape(j, 0) * vk[0];
1811  if (vshape.Width() == 3)
1812  {
1813  I(k, j) += vshape(j, 1) * t3(1);
1814  I(k, j) += vshape(j, 2) * t3(2);
1815  }
1816  }
1817  }
1818  }
1819 }
1820 
1821 const double ND_R2D_SegmentElement::tk[4] = { 1.,0., 0.,1. };
1822 
1824  const int cb_type,
1825  const int ob_type)
1826  : VectorFiniteElement(1, Geometry::SEGMENT, 2 * p + 1, p,
1827  H_CURL, FunctionSpace::Pk),
1828  dof2tk(dof),
1829  cbasis1d(poly1d.GetBasis(p, VerifyClosed(cb_type))),
1830  obasis1d(poly1d.GetBasis(p - 1, VerifyOpen(ob_type)))
1831 {
1832  // Override default dimensions for VectorFiniteElements
1833  vdim = 2;
1834  cdim = 1;
1835 
1836  const double *cp = poly1d.ClosedPoints(p, cb_type);
1837  const double *op = poly1d.OpenPoints(p - 1, ob_type);
1838 
1839 #ifndef MFEM_THREAD_SAFE
1840  shape_cx.SetSize(p + 1);
1841  shape_ox.SetSize(p);
1842  dshape_cx.SetSize(p + 1);
1843 #endif
1844 
1845  dof_map.SetSize(dof);
1846 
1847  int o = 0;
1848  // nodes
1849  // (0)
1850  Nodes.IntPoint(o).x = cp[0]; // z-directed
1851  dof_map[p] = o; dof2tk[o++] = 1;
1852 
1853  // (1)
1854  Nodes.IntPoint(o).x = cp[p]; // z-directed
1855  dof_map[2*p] = o; dof2tk[o++] = 1;
1856 
1857  // interior
1858  // x-components
1859  for (int i = 0; i < p; i++)
1860  {
1861  Nodes.IntPoint(o).x = op[i];
1862  dof_map[i] = o; dof2tk[o++] = 0;
1863  }
1864  // z-components
1865  for (int i = 1; i < p; i++)
1866  {
1867  Nodes.IntPoint(o).x = cp[i];
1868  dof_map[p+i] = o; dof2tk[o++] = 1;
1869  }
1870 }
1871 
1873  DenseMatrix &shape) const
1874 {
1875  const int p = order;
1876 
1877 #ifdef MFEM_THREAD_SAFE
1878  Vector shape_cx(p + 1), shape_ox(p);
1879 #endif
1880 
1881  cbasis1d.Eval(ip.x, shape_cx);
1882  obasis1d.Eval(ip.x, shape_ox);
1883 
1884  int o = 0;
1885  // x-components
1886  for (int i = 0; i < p; i++)
1887  {
1888  int idx = dof_map[o++];
1889  shape(idx,0) = shape_ox(i);
1890  shape(idx,1) = 0.;
1891  }
1892  // z-components
1893  for (int i = 0; i <= p; i++)
1894  {
1895  int idx = dof_map[o++];
1896  shape(idx,0) = 0.;
1897  shape(idx,1) = shape_cx(i);
1898  }
1899 }
1900 
1902  DenseMatrix &shape) const
1903 {
1904  CalcVShape(Trans.GetIntPoint(), shape);
1905  const DenseMatrix & JI = Trans.InverseJacobian();
1906  MFEM_ASSERT(JI.Width() == 1 && JI.Height() == 1,
1907  "ND_R2D_SegmentElement cannot be embedded in "
1908  "2 or 3 dimensional spaces");
1909  for (int i=0; i<dof; i++)
1910  {
1911  shape(i, 0) *= JI(0,0);
1912  }
1913 }
1914 
1916  DenseMatrix &curl_shape) const
1917 {
1918  const int p = order;
1919 
1920 #ifdef MFEM_THREAD_SAFE
1921  Vector shape_cx(p + 1), shape_ox(p);
1922  Vector dshape_cx(p + 1);
1923 #endif
1924 
1925  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
1926  obasis1d.Eval(ip.x, shape_ox);
1927 
1928  int o = 0;
1929  // x-components
1930  for (int i = 0; i < p; i++)
1931  {
1932  int idx = dof_map[o++];
1933  curl_shape(idx,0) = 0.;
1934  }
1935  // z-components
1936  for (int i = 0; i <= p; i++)
1937  {
1938  int idx = dof_map[o++];
1939  curl_shape(idx,0) = -dshape_cx(i);
1940  }
1941 }
1942 
1943 void ND_R2D_SegmentElement::LocalInterpolation(const VectorFiniteElement &cfe,
1945  DenseMatrix &I) const
1946 {
1947  double vk[Geometry::MaxDim]; vk[1] = 0.0; vk[2] = 0.0;
1948  Vector xk(vk, dim);
1949  IntegrationPoint ip;
1950  DenseMatrix vshape(cfe.GetDof(), vdim);
1951 
1952  double * tk_ptr = const_cast<double*>(tk);
1953 
1954  I.SetSize(dof, vshape.Height());
1955 
1956  // assuming Trans is linear; this should be ok for all refinement types
1957  Trans.SetIntPoint(&Geometries.GetCenter(geom_type));
1958  const DenseMatrix &J = Trans.Jacobian();
1959  for (int k = 0; k < dof; k++)
1960  {
1961  Vector t1(&tk_ptr[dof2tk[k] * 2], 1);
1962  Vector t2(&tk_ptr[dof2tk[k] * 2], 2);
1963 
1964  Trans.Transform(Nodes.IntPoint(k), xk);
1965  ip.Set3(vk);
1966  cfe.CalcVShape(ip, vshape);
1967  // xk = J t_k
1968  J.Mult(t1, vk);
1969  // I_k = vshape_k.J.t_k, k=1,...,Dof
1970  for (int j = 0; j < vshape.Height(); j++)
1971  {
1972  double Ikj = 0.;
1973  for (int i = 0; i < dim; i++)
1974  {
1975  Ikj += vshape(j, i) * vk[i];
1976  }
1977  Ikj += vshape(j, 1) * t2(1);
1978  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1979  }
1980  }
1981 }
1982 
1985  Vector &dofs) const
1986 {
1987  double data[3];
1988  Vector vk1(data, 1);
1989  Vector vk2(data, 2);
1990  Vector vk3(data, 3);
1991 
1992  double * tk_ptr = const_cast<double*>(tk);
1993 
1994  for (int k = 0; k < dof; k++)
1995  {
1996  Trans.SetIntPoint(&Nodes.IntPoint(k));
1997 
1998  vc.Eval(vk3, Trans, Nodes.IntPoint(k));
1999  // dof_k = vk^t J tk
2000  Vector t1(&tk_ptr[dof2tk[k] * 2], 1);
2001  Vector t2(&tk_ptr[dof2tk[k] * 2], 2);
2002 
2003  dofs(k) = Trans.Jacobian().InnerProduct(t1, vk2) + t2(1) * vk3(2);
2004  }
2005 
2006 }
2007 
2009  const double *tk_fe)
2010  : VectorFiniteElement(2, G, Do, p,
2011  H_CURL, FunctionSpace::Pk),
2012  tk(tk_fe),
2013  dof_map(dof),
2014  dof2tk(dof)
2015 {
2016  // Override default types for VectorFiniteElements
2017  deriv_type = CURL;
2020 
2021  // Override default dimensions for VectorFiniteElements
2022  vdim = 3;
2023  cdim = 3;
2024 }
2025 
2027  DenseMatrix &shape) const
2028 {
2029  CalcVShape(Trans.GetIntPoint(), shape);
2030  const DenseMatrix & JI = Trans.InverseJacobian();
2031  MFEM_ASSERT(JI.Width() == 2 && JI.Height() == 2,
2032  "ND_R2D_FiniteElement cannot be embedded in "
2033  "3 dimensional spaces");
2034  for (int i=0; i<dof; i++)
2035  {
2036  double sx = shape(i, 0);
2037  double sy = shape(i, 1);
2038  shape(i, 0) = sx * JI(0, 0) + sy * JI(1, 0);
2039  shape(i, 1) = sx * JI(0, 1) + sy * JI(1, 1);
2040  }
2041 }
2042 
2044  DenseMatrix &curl_shape) const
2045 {
2046  CalcCurlShape(Trans.GetIntPoint(), curl_shape);
2047  const DenseMatrix & J = Trans.Jacobian();
2048  MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
2049  "ND_R2D_FiniteElement cannot be embedded in "
2050  "3 dimensional spaces");
2051  for (int i=0; i<dof; i++)
2052  {
2053  double sx = curl_shape(i, 0);
2054  double sy = curl_shape(i, 1);
2055  curl_shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
2056  curl_shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
2057  }
2058  curl_shape *= (1.0 / Trans.Weight());
2059 }
2060 
2061 void ND_R2D_FiniteElement::LocalInterpolation(
2062  const VectorFiniteElement &cfe,
2064  DenseMatrix &I) const
2065 {
2066  double vk[Geometry::MaxDim]; vk[2] = 0.0;
2067  Vector xk(vk, dim);
2068  IntegrationPoint ip;
2069 #ifdef MFEM_THREAD_SAFE
2070  DenseMatrix vshape(cfe.GetDof(), vdim);
2071 #else
2072  vshape.SetSize(cfe.GetDof(), vdim);
2073 #endif
2074 
2075  double * tk_ptr = const_cast<double*>(tk);
2076 
2077  I.SetSize(dof, vshape.Height());
2078 
2079  // assuming Trans is linear; this should be ok for all refinement types
2080  Trans.SetIntPoint(&Geometries.GetCenter(geom_type));
2081  const DenseMatrix &J = Trans.Jacobian();
2082  for (int k = 0; k < dof; k++)
2083  {
2084  Vector t2(&tk_ptr[dof2tk[k] * 3], 2);
2085  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
2086 
2087  Trans.Transform(Nodes.IntPoint(k), xk);
2088  ip.Set3(vk);
2089  cfe.CalcVShape(ip, vshape);
2090  // xk = J t_k
2091  J.Mult(t2, vk);
2092  // I_k = vshape_k.J.t_k, k=1,...,Dof
2093  for (int j = 0; j < vshape.Height(); j++)
2094  {
2095  double Ikj = 0.;
2096  for (int i = 0; i < dim; i++)
2097  {
2098  Ikj += vshape(j, i) * vk[i];
2099  }
2100  Ikj += vshape(j, 2) * t3(2);
2101  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
2102  }
2103  }
2104 }
2105 
2107  DenseMatrix &R) const
2108 {
2109  double pt_data[Geometry::MaxDim];
2110  IntegrationPoint ip;
2111  Vector pt(pt_data, dim);
2112 
2113 #ifdef MFEM_THREAD_SAFE
2115 #endif
2116 
2117  double * tk_ptr = const_cast<double*>(tk);
2118 
2119  Trans.SetIntPoint(&Geometries.GetCenter(geom_type));
2120  const DenseMatrix &Jinv = Trans.InverseJacobian();
2121  for (int j = 0; j < dof; j++)
2122  {
2123  Vector t2(&tk_ptr[dof2tk[j] * 3], 2);
2124  Vector t3(&tk_ptr[dof2tk[j] * 3], 3);
2125 
2127  ip.Set(pt_data, dim);
2128  if (Geometries.CheckPoint(geom_type, ip)) // do we need an epsilon here?
2129  {
2130  CalcVShape(ip, vshape);
2131  Jinv.Mult(t2, pt_data);
2132  for (int k = 0; k < dof; k++)
2133  {
2134  double R_jk = 0.0;
2135  for (int d = 0; d < dim; d++)
2136  {
2137  R_jk += vshape(k,d)*pt_data[d];
2138  }
2139  R_jk += vshape(k, 2) * t3(2);
2140  R(j,k) = R_jk;
2141  }
2142  }
2143  else
2144  {
2145  // Set the whole row to avoid valgrind warnings in R.Threshold().
2146  R.SetRow(j, infinity());
2147  }
2148  }
2149  R.Threshold(1e-12);
2150 }
2151 
2154  Vector &dofs) const
2155 {
2156  double data[3];
2157  Vector vk2(data, 2);
2158  Vector vk3(data, 3);
2159 
2160  double * tk_ptr = const_cast<double*>(tk);
2161 
2162  for (int k = 0; k < dof; k++)
2163  {
2164  Trans.SetIntPoint(&Nodes.IntPoint(k));
2165 
2166  vc.Eval(vk3, Trans, Nodes.IntPoint(k));
2167  // dof_k = vk^t J tk
2168  Vector t2(&tk_ptr[dof2tk[k] * 3], 2);
2169  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
2170 
2171  dofs(k) = Trans.Jacobian().InnerProduct(t2, vk2) + t3(2) * vk3(2);
2172  }
2173 
2174 }
2175 
2178  DenseMatrix &I) const
2179 {
2180  if (fe.GetRangeType() == SCALAR)
2181  {
2182  double vk[Geometry::MaxDim];
2183  Vector shape(fe.GetDof());
2184 
2185  double * tk_ptr = const_cast<double*>(tk);
2186 
2187  I.SetSize(dof, vdim*fe.GetDof());
2188  for (int k = 0; k < dof; k++)
2189  {
2190  const IntegrationPoint &ip = Nodes.IntPoint(k);
2191 
2192  Vector t2(&tk_ptr[dof2tk[k] * 3], 2);
2193  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
2194 
2195  fe.CalcShape(ip, shape);
2196  Trans.SetIntPoint(&ip);
2197  // Transform ND edge tengents from reference to physical space
2198  // vk = J tk
2199  Trans.Jacobian().Mult(t2, vk);
2200  vk[2] = t3[2];
2201  if (fe.GetMapType() == INTEGRAL)
2202  {
2203  double w = 1.0/Trans.Weight();
2204  for (int d = 0; d < vdim; d++)
2205  {
2206  vk[d] *= w;
2207  }
2208  }
2209 
2210  for (int j = 0; j < shape.Size(); j++)
2211  {
2212  double s = shape(j);
2213  if (fabs(s) < 1e-12)
2214  {
2215  s = 0.0;
2216  }
2217  // Project scalar basis function multiplied by each coordinate
2218  // direction onto the transformed edge tangents
2219  for (int d = 0; d < vdim; d++)
2220  {
2221  I(k, j + d*shape.Size()) = s*vk[d];
2222  }
2223  }
2224  }
2225  }
2226  else
2227  {
2228  double vk[Geometry::MaxDim];
2229  DenseMatrix vshape(fe.GetDof(), fe.GetVDim());
2230 
2231  double * tk_ptr = const_cast<double*>(tk);
2232 
2233  I.SetSize(dof, fe.GetDof());
2234  for (int k = 0; k < dof; k++)
2235  {
2236  const IntegrationPoint &ip = Nodes.IntPoint(k);
2237 
2238  Vector t2(&tk_ptr[dof2tk[k] * 3], 2);
2239  Vector t3(&tk_ptr[dof2tk[k] * 3], 3);
2240 
2241  Trans.SetIntPoint(&ip);
2242  // Transform ND edge tangents from reference to physical space
2243  // vk = J tk
2244  Trans.Jacobian().Mult(t2, vk);
2245  // Compute fe basis functions in physical space
2246  fe.CalcVShape(Trans, vshape);
2247  // Project fe basis functions onto transformed edge tangents
2248  for (int j=0; j<vshape.Height(); j++)
2249  {
2250  I(k, j) = 0.0;
2251  for (int i=0; i<2; i++)
2252  {
2253  I(k, j) += vshape(j, i) * vk[i];
2254  }
2255  if (vshape.Width() == 3)
2256  {
2257  I(k, j) += vshape(j, 2) * t3(2);
2258  }
2259  }
2260  }
2261  }
2262 }
2263 
2266  DenseMatrix &grad) const
2267 {
2268  MFEM_ASSERT(fe.GetMapType() == VALUE, "");
2269 
2270  DenseMatrix dshape(fe.GetDof(), fe.GetDim());
2271  Vector grad_k(fe.GetDof());
2272 
2273  grad.SetSize(dof, fe.GetDof());
2274  for (int k = 0; k < dof; k++)
2275  {
2276  fe.CalcDShape(Nodes.IntPoint(k), dshape);
2277  dshape.Mult(tk + dof2tk[k]*vdim, grad_k);
2278  for (int j = 0; j < grad_k.Size(); j++)
2279  {
2280  grad(k,j) = (fabs(grad_k(j)) < 1e-12) ? 0.0 : grad_k(j);
2281  }
2282  }
2283 }
2284 
2285 const double ND_R2D_TriangleElement::tk_t[15] =
2286 { 1.,0.,0., -1.,1.,0., 0.,-1.,0., 0.,1.,0., 0.,0.,1. };
2287 
2289  const int cb_type)
2290  : ND_R2D_FiniteElement(p, Geometry::TRIANGLE, ((3*p + 1)*(p + 2))/2, tk_t),
2291  ND_FE(p),
2292  H1_FE(p, cb_type)
2293 {
2294  int pm1 = p - 1, pm2 = p - 2;
2295 
2296 #ifndef MFEM_THREAD_SAFE
2297  nd_shape.SetSize(ND_FE.GetDof(), 2);
2298  h1_shape.SetSize(H1_FE.GetDof());
2299  nd_dshape.SetSize(ND_FE.GetDof(), 1);
2300  h1_dshape.SetSize(H1_FE.GetDof(), 2);
2301 #endif
2302 
2303  int o = 0;
2304  int n = 0;
2305  int h = 0;
2306  // Three nodes
2307  dof_map[o] = -1 - h++; dof2tk[o++] = 4;
2308  dof_map[o] = -1 - h++; dof2tk[o++] = 4;
2309  dof_map[o] = -1 - h++; dof2tk[o++] = 4;
2310 
2311  // Three edges
2312  for (int e=0; e<3; e++)
2313  {
2314  // Dofs in the plane
2315  for (int i=0; i<p; i++)
2316  {
2317  dof_map[o] = n++; dof2tk[o++] = e;
2318  }
2319  // z-directed dofs
2320  for (int i=0; i<pm1; i++)
2321  {
2322  dof_map[o] = -1 - h++; dof2tk[o++] = 4;
2323  }
2324  }
2325 
2326  // Interior dofs in the plane
2327  for (int j = 0; j <= pm2; j++)
2328  for (int i = 0; i + j <= pm2; i++)
2329  {
2330  dof_map[o] = n++; dof2tk[o++] = 0;
2331  dof_map[o] = n++; dof2tk[o++] = 3;
2332  }
2333 
2334  // Interior z-directed dofs
2335  for (int j = 0; j < pm1; j++)
2336  for (int i = 0; i + j < pm2; i++)
2337  {
2338  dof_map[o] = -1 - h++; dof2tk[o++] = 4;
2339  }
2340 
2341  MFEM_VERIFY(n == ND_FE.GetDof(),
2342  "ND_R2D_Triangle incorrect number of ND dofs.");
2343  MFEM_VERIFY(h == H1_FE.GetDof(),
2344  "ND_R2D_Triangle incorrect number of H1 dofs.");
2345  MFEM_VERIFY(o == GetDof(),
2346  "ND_R2D_Triangle incorrect number of dofs.");
2347 
2348  const IntegrationRule & nd_Nodes = ND_FE.GetNodes();
2349  const IntegrationRule & h1_Nodes = H1_FE.GetNodes();
2350 
2351  for (int i=0; i<dof; i++)
2352  {
2353  int idx = dof_map[i];
2354  if (idx >= 0)
2355  {
2356  const IntegrationPoint & ip = nd_Nodes.IntPoint(idx);
2357  Nodes.IntPoint(i).Set2(ip.x, ip.y);
2358  }
2359  else
2360  {
2361  const IntegrationPoint & ip = h1_Nodes.IntPoint(-idx-1);
2362  Nodes.IntPoint(i).Set2(ip.x, ip.y);
2363  }
2364  }
2365 }
2366 
2368  DenseMatrix &shape) const
2369 {
2370 #ifdef MFEM_THREAD_SAFE
2371  DenseMatrix nd_shape(ND_FE.GetDof(), 2);
2372  Vector h1_shape(H1_FE.GetDof());
2373 #endif
2374 
2375  ND_FE.CalcVShape(ip, nd_shape);
2376  H1_FE.CalcShape(ip, h1_shape);
2377 
2378  for (int i=0; i<dof; i++)
2379  {
2380  int idx = dof_map[i];
2381  if (idx >= 0)
2382  {
2383  shape(i, 0) = nd_shape(idx, 0);
2384  shape(i, 1) = nd_shape(idx, 1);
2385  shape(i, 2) = 0.0;
2386  }
2387  else
2388  {
2389  shape(i, 0) = 0.0;
2390  shape(i, 1) = 0.0;
2391  shape(i, 2) = h1_shape(-idx-1);
2392  }
2393  }
2394 }
2395 
2397  DenseMatrix &curl_shape) const
2398 {
2399 #ifdef MFEM_THREAD_SAFE
2400  DenseMatrix nd_dshape(ND_FE.GetDof(), 1);
2401  DenseMatrix h1_dshape(H1_FE.GetDof(), 2);
2402 #endif
2403 
2404  ND_FE.CalcCurlShape(ip, nd_dshape);
2405  H1_FE.CalcDShape(ip, h1_dshape);
2406 
2407  for (int i=0; i<dof; i++)
2408  {
2409  int idx = dof_map[i];
2410  if (idx >= 0)
2411  {
2412  curl_shape(i, 0) = 0.0;
2413  curl_shape(i, 1) = 0.0;
2414  curl_shape(i, 2) = nd_dshape(idx, 0);
2415  }
2416  else
2417  {
2418  curl_shape(i, 0) = h1_dshape(-idx-1, 1);
2419  curl_shape(i, 1) = -h1_dshape(-idx-1, 0);
2420  curl_shape(i, 2) = 0.0;
2421  }
2422  }
2423 }
2424 
2425 
2426 const double ND_R2D_QuadrilateralElement::tk_q[15] =
2427 { 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,-1.,0., 0.,0.,1. };
2428 
2430  const int cb_type,
2431  const int ob_type)
2432  : ND_R2D_FiniteElement(p, Geometry::SQUARE, ((3*p + 1)*(p + 1)), tk_q),
2433  cbasis1d(poly1d.GetBasis(p, VerifyClosed(cb_type))),
2434  obasis1d(poly1d.GetBasis(p - 1, VerifyOpen(ob_type)))
2435 {
2436  const double *cp = poly1d.ClosedPoints(p, cb_type);
2437  const double *op = poly1d.OpenPoints(p - 1, ob_type);
2438  const int dofx = p*(p+1);
2439  const int dofy = p*(p+1);
2440  const int dofxy = dofx+dofy;
2441 
2442 #ifndef MFEM_THREAD_SAFE
2443  shape_cx.SetSize(p + 1);
2444  shape_ox.SetSize(p);
2445  shape_cy.SetSize(p + 1);
2446  shape_oy.SetSize(p);
2447  dshape_cx.SetSize(p + 1);
2448  dshape_cy.SetSize(p + 1);
2449 #endif
2450 
2451  dof_map.SetSize(dof);
2452 
2453  int o = 0;
2454  // nodes
2455  dof_map[dofxy] = o++; // (0)
2456  dof_map[dofxy+p] = o++; // (1)
2457  dof_map[dof-1] = o++; // (2)
2458  dof_map[dof-p-1] = o++; // (3)
2459 
2460  // edges
2461  for (int i = 0; i < p; i++) // (0,1) - x-directed
2462  {
2463  dof_map[i + 0*p] = o++;
2464  }
2465  for (int i = 1; i < p; i++) // (0,1) - z-directed
2466  {
2467  dof_map[dofxy + i + 0*(p+1)] = o++;
2468  }
2469  for (int j = 0; j < p; j++) // (1,2) - y-directed
2470  {
2471  dof_map[dofx + p + j*(p + 1)] = o++;
2472  }
2473  for (int j = 1; j < p; j++) // (1,2) - z-directed
2474  {
2475  dof_map[dofxy + p + j*(p + 1)] = o++;
2476  }
2477  for (int i = 0; i < p; i++) // (2,3) - x-directed
2478  {
2479  dof_map[(p - 1 - i) + p*p] = -1 - (o++);
2480  }
2481  for (int i = 1; i < p; i++) // (2,3) - z-directed
2482  {
2483  dof_map[dofxy + (p - i) + p*(p + 1)] = o++;
2484  }
2485  for (int j = 0; j < p; j++) // (3,0) - y-directed
2486  {
2487  dof_map[dofx + 0 + (p - 1 - j)*(p + 1)] = -1 - (o++);
2488  }
2489  for (int j = 1; j < p; j++) // (3,0) - z-directed
2490  {
2491  dof_map[dofxy + (p - j)*(p + 1)] = o++;
2492  }
2493 
2494  // interior
2495  // x-components
2496  for (int j = 1; j < p; j++)
2497  for (int i = 0; i < p; i++)
2498  {
2499  dof_map[i + j*p] = o++;
2500  }
2501  // y-components
2502  for (int j = 0; j < p; j++)
2503  for (int i = 1; i < p; i++)
2504  {
2505  dof_map[dofx + i + j*(p + 1)] = o++;
2506  }
2507  // z-components
2508  for (int j = 1; j < p; j++)
2509  for (int i = 1; i < p; i++)
2510  {
2511  dof_map[dofxy + i + j*(p + 1)] = o++;
2512  }
2513 
2514  // set dof2tk and Nodes
2515  o = 0;
2516  // x-components
2517  for (int j = 0; j <= p; j++)
2518  for (int i = 0; i < p; i++)
2519  {
2520  int idx;
2521  if ((idx = dof_map[o++]) < 0)
2522  {
2523  dof2tk[idx = -1 - idx] = 2;
2524  }
2525  else
2526  {
2527  dof2tk[idx] = 0;
2528  }
2529  Nodes.IntPoint(idx).Set2(op[i], cp[j]);
2530  }
2531  // y-components
2532  for (int j = 0; j < p; j++)
2533  for (int i = 0; i <= p; i++)
2534  {
2535  int idx;
2536  if ((idx = dof_map[o++]) < 0)
2537  {
2538  dof2tk[idx = -1 - idx] = 3;
2539  }
2540  else
2541  {
2542  dof2tk[idx] = 1;
2543  }
2544  Nodes.IntPoint(idx).Set2(cp[i], op[j]);
2545  }
2546  // z-components
2547  for (int j = 0; j <= p; j++)
2548  for (int i = 0; i <= p; i++)
2549  {
2550  int idx = dof_map[o++];
2551  dof2tk[idx] = 4;
2552  Nodes.IntPoint(idx).Set2(cp[i], cp[j]);
2553  }
2554 }
2555 
2557  DenseMatrix &shape) const
2558 {
2559  const int p = order;
2560 
2561 #ifdef MFEM_THREAD_SAFE
2562  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
2563 #endif
2564 
2565  cbasis1d.Eval(ip.x, shape_cx);
2566  obasis1d.Eval(ip.x, shape_ox);
2567  cbasis1d.Eval(ip.y, shape_cy);
2568  obasis1d.Eval(ip.y, shape_oy);
2569 
2570  int o = 0;
2571  // x-components
2572  for (int j = 0; j <= p; j++)
2573  for (int i = 0; i < p; i++)
2574  {
2575  int idx, s;
2576  if ((idx = dof_map[o++]) < 0)
2577  {
2578  idx = -1 - idx, s = -1;
2579  }
2580  else
2581  {
2582  s = +1;
2583  }
2584  shape(idx,0) = s*shape_ox(i)*shape_cy(j);
2585  shape(idx,1) = 0.;
2586  shape(idx,2) = 0.;
2587  }
2588  // y-components
2589  for (int j = 0; j < p; j++)
2590  for (int i = 0; i <= p; i++)
2591  {
2592  int idx, s;
2593  if ((idx = dof_map[o++]) < 0)
2594  {
2595  idx = -1 - idx, s = -1;
2596  }
2597  else
2598  {
2599  s = +1;
2600  }
2601  shape(idx,0) = 0.;
2602  shape(idx,1) = s*shape_cx(i)*shape_oy(j);
2603  shape(idx,2) = 0.;
2604  }
2605  // z-components
2606  for (int j = 0; j <= p; j++)
2607  for (int i = 0; i <= p; i++)
2608  {
2609  int idx = dof_map[o++];
2610  shape(idx,0) = 0.;
2611  shape(idx,1) = 0.;
2612  shape(idx,2) = shape_cx(i)*shape_cy(j);
2613  }
2614 }
2615 
2617  DenseMatrix &curl_shape) const
2618 {
2619  const int p = order;
2620 
2621 #ifdef MFEM_THREAD_SAFE
2622  Vector shape_cx(p + 1), shape_ox(p), shape_cy(p + 1), shape_oy(p);
2623  Vector dshape_cx(p + 1), dshape_cy(p + 1);
2624 #endif
2625 
2626  cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
2627  obasis1d.Eval(ip.x, shape_ox);
2628  cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
2629  obasis1d.Eval(ip.y, shape_oy);
2630 
2631  int o = 0;
2632  // x-components
2633  for (int j = 0; j <= p; j++)
2634  for (int i = 0; i < p; i++)
2635  {
2636  int idx, s;
2637  if ((idx = dof_map[o++]) < 0)
2638  {
2639  idx = -1 - idx, s = -1;
2640  }
2641  else
2642  {
2643  s = +1;
2644  }
2645  curl_shape(idx,0) = 0.;
2646  curl_shape(idx,1) = 0.;
2647  curl_shape(idx,2) = -s*shape_ox(i)*dshape_cy(j);
2648  }
2649  // y-components
2650  for (int j = 0; j < p; j++)
2651  for (int i = 0; i <= p; i++)
2652  {
2653  int idx, s;
2654  if ((idx = dof_map[o++]) < 0)
2655  {
2656  idx = -1 - idx, s = -1;
2657  }
2658  else
2659  {
2660  s = +1;
2661  }
2662  curl_shape(idx,0) = 0.;
2663  curl_shape(idx,1) = 0.;
2664  curl_shape(idx,2) = s*dshape_cx(i)*shape_oy(j);
2665  }
2666  // z-components
2667  for (int j = 0; j <= p; j++)
2668  for (int i = 0; i <= p; i++)
2669  {
2670  int idx = dof_map[o++];
2671  curl_shape(idx,0) = shape_cx(i)*dshape_cy(j);
2672  curl_shape(idx,1) = -dshape_cx(i)*shape_cy(j);
2673  curl_shape(idx,2) = 0.;
2674  }
2675 }
2676 
2677 }
Abstract class for all finite elements.
Definition: fe_base.hpp:232
void Set(const double *p, const int dim)
Definition: intrules.hpp:37
ND_SegmentElement(const int p, const int ob_type=BasisType::GaussLegendre)
Construct the ND_SegmentElement of order p and open BasisType ob_type.
Definition: fe_nd.cpp:1202
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_nd.cpp:1528
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_nd.cpp:1479
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
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
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_nd.cpp:1714
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
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_nd.cpp:654
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 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_nd.cpp:1219
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
ND_R1D_PointElement(int p)
Construct the ND_R1D_PointElement.
Definition: fe_nd.cpp:1517
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 ...
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:3963
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:512
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe_nd.cpp:244
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_nd.cpp:717
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_nd.cpp:1872
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
ND_R2D_TriangleElement(const int p, const int cb_type=BasisType::GaussLobatto)
Construct the ND_R2D_TriangleElement of order p.
Definition: fe_nd.cpp:2288
ND_R2D_SegmentElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the ND_R2D_SegmentElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_nd.cpp:1823
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 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_nd.cpp:2026
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
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe_nd.cpp:578
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(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_nd.cpp:1123
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
ND_WedgeElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_nd.cpp:1230
Geometry Geometries
Definition: fe.cpp:49
ND_TriangleElement(const int p)
Construct the ND_TriangleElement of order p.
Definition: fe_nd.cpp:1042
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
int cdim
Dimension of curl for vector-valued basis functions.
Definition: fe_base.hpp:237
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 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_nd.cpp:1610
ND_R2D_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the ND_R2D_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_typ...
Definition: fe_nd.cpp:2429
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_nd.cpp:1157
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
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_nd.cpp:1445
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
Definition: fe_nd.cpp:2043
ND_HexahedronElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the ND_HexahedronElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_nd.cpp:25
Poly_1D::Basis & basis1d
Definition: fe_base.hpp:1166
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:555
ND_R2D_FiniteElement(int p, Geometry::Type G, int Do, const double *tk_fe)
Definition: fe_nd.cpp:2008
ND_TetrahedronElement(const int p)
Construct the ND_TetrahedronElement of order p.
Definition: fe_nd.cpp:780
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_nd.cpp:2616
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 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_nd.cpp:970
const double * tk
Definition: fe_nd.hpp:552
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_nd.cpp:928
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
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_nd.cpp:2152
ND_R1D_SegmentElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the ND_R1D_SegmentElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_nd.cpp:1546
void Set2(const double x1, const double x2)
Definition: intrules.hpp:80
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_nd.cpp:1663
int GetVDim()
Returns dimension of the vector.
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:173
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_nd.cpp:2556
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 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_nd.cpp:399
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
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_nd.cpp:309
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:70
IntegrationRule Nodes
Definition: fe_base.hpp:245
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_nd.cpp:1915
ND_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the ND_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type...
Definition: fe_nd.cpp:487
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:310
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:684
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
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:533
Class for integration point with weight.
Definition: intrules.hpp:25
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 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
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_nd.cpp:2396
RefCoord t[3]
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_nd.cpp:2106
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 Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Definition: fe_nd.cpp:1983
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
RefCoord s[3]
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_nd.cpp:2367
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Array< int > dof_map
Definition: fe_nd.hpp:553
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:388
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
Definition: fe_nd.cpp:1703
int order
Order/degree of the shape functions.
Definition: fe_base.hpp:242
Implements CalcCurlShape methods.
Definition: fe_base.hpp:296
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient 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_nd.cpp:2264