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