libMesh
fe_monomial_vec.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 // Local includes
19 #include "libmesh/dof_map.h"
20 #include "libmesh/elem.h"
21 #include "libmesh/fe.h"
22 #include "libmesh/fe_interface.h"
23 #include "libmesh/fe_macro.h"
24 #include "libmesh/tensor_value.h"
25 
26 
27 namespace libMesh
28 {
29 
30 
35 
36 
37 // ------------------------------------------------------------
38 // Vector monomial specific implementations
39 
40 // Anonymous namespace for local helper functions
41 namespace
42 {
43 void
44 monomial_vec_nodal_soln(const Elem * elem,
45  const Order order,
46  const std::vector<Number> & elem_soln,
47  const int dim,
48  std::vector<Number> & nodal_soln,
49  const bool add_p_level)
50 {
51  const unsigned int n_nodes = elem->n_nodes();
52 
53  const ElemType elem_type = elem->type();
54 
55  nodal_soln.resize(dim * n_nodes);
56 
57  const Order totalorder = order + add_p_level*elem->p_level();
58 
59  switch (totalorder)
60  {
61  // Constant shape functions
62  case CONSTANT:
63  {
64  libmesh_assert_equal_to(elem_soln.size(), static_cast<unsigned int>(dim));
65  switch (dim)
66  {
67  case 2:
68  case 3:
69  {
70  for (unsigned int n = 0; n < n_nodes; n++)
71  std::copy(elem_soln.begin(), elem_soln.end(), nodal_soln.begin() + dim*n);
72  return;
73  }
74  default:
75  libmesh_error_msg(
76  "The monomial_vec_nodal_soln helper should only be called for 2 and 3 dimensions");
77  }
78  }
79 
80  // For other orders, do interpolation at the nodes
81  // explicitly.
82  default:
83  {
84  // FEType object to be passed to various FEInterface functions below.
85  FEType fe_type(order, MONOMIAL);
86 
87  const unsigned int n_sf = FEInterface::n_shape_functions(fe_type, elem);
88 
89  std::vector<Point> refspace_nodes;
90  FEBase::get_refspace_nodes(elem_type, refspace_nodes);
91  libmesh_assert_equal_to(refspace_nodes.size(), n_nodes);
92  libmesh_assert_equal_to(elem_soln.size(), n_sf * dim);
93 
94  // Zero before summation
95  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
96 
97  for (unsigned int d = 0; d < static_cast<unsigned int>(dim); d++)
98  for (unsigned int n = 0; n < n_nodes; n++)
99  // u_i = Sum (alpha_i phi_i)
100  for (unsigned int i = 0; i < n_sf; i++)
101  nodal_soln[d + dim * n] += elem_soln[d + dim * i] *
102  FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
103 
104  return;
105  } // default
106  } // switch
107 }
108 } // anonymous namespace
109 
110 // Do full-specialization for every dimension, instead
111 // of explicit instantiation at the end of this file.
112 // This could be macro-ified so that it fits on one line...
114 LIBMESH_FE_NODAL_SOLN_DIM(MONOMIAL_VEC, (FE<1, MONOMIAL>::nodal_soln), 1)
115 
116 template <>
117 void
118 FE<2, MONOMIAL_VEC>::nodal_soln(const Elem * elem,
119  const Order order,
120  const std::vector<Number> & elem_soln,
121  std::vector<Number> & nodal_soln,
122  const bool add_p_level,
123  const unsigned)
124 {
125  monomial_vec_nodal_soln(elem, order, elem_soln, 2 /*dim*/, nodal_soln, add_p_level);
126 }
127 
128 template <>
129 void
131  const Order order,
132  const std::vector<Number> & elem_soln,
133  std::vector<Number> & nodal_soln,
134  const bool add_p_level,
135  const unsigned)
136 {
137  monomial_vec_nodal_soln(elem, order, elem_soln, 3 /*dim*/, nodal_soln, add_p_level);
138 }
139 
141 
142 
143 // Specialize for shape function routines by leveraging scalar MONOMIAL elements
144 
145 // 0-D
146 template <>
149  const Order order,
150  const unsigned int i,
151  const Point & p)
152 {
153  Real value = FE<0, MONOMIAL>::shape(type, order, i, p);
155 }
156 template <>
159  const Order order,
160  const unsigned int i,
161  const unsigned int j,
162  const Point & p)
163 {
164  Real value = FE<0, MONOMIAL>::shape_deriv(type, order, i, j, p);
166 }
167 
168 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
169 
170 template <>
173  const Order order,
174  const unsigned int i,
175  const unsigned int j,
176  const Point & p)
177 {
178  Real value = FE<0, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
180 }
181 #endif
182 
183 // 1-D
184 template <>
187  const Order order,
188  const unsigned int i,
189  const Point & p)
190 {
191  Real value = FE<1, MONOMIAL>::shape(type, order, i, p);
193 }
194 template <>
197  const Order order,
198  const unsigned int i,
199  const unsigned int j,
200  const Point & p)
201 {
202  Real value = FE<1, MONOMIAL>::shape_deriv(type, order, i, j, p);
204 }
205 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
206 template <>
209  const Order order,
210  const unsigned int i,
211  const unsigned int j,
212  const Point & p)
213 {
214  Real value = FE<1, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
216 }
217 
218 #endif
219 
220 // 2-D
221 template <>
224  const Order order,
225  const unsigned int i,
226  const Point & p)
227 {
228  Real value = FE<2, MONOMIAL>::shape(type, order, i / 2, p);
229 
230  switch (i % 2)
231  {
232  case 0:
234 
235  case 1:
236  return libMesh::RealVectorValue(Real(0), value);
237 
238  default:
239  libmesh_error_msg("i%2 must be either 0 or 1!");
240  }
241 
242  // dummy
243  return libMesh::RealVectorValue();
244 }
245 template <>
248  const Order order,
249  const unsigned int i,
250  const unsigned int j,
251  const Point & p)
252 {
253  Real value = FE<2, MONOMIAL>::shape_deriv(type, order, i / 2, j, p);
254 
255  switch (i % 2)
256  {
257  case 0:
259 
260  case 1:
261  return libMesh::RealVectorValue(Real(0), value);
262 
263  default:
264  libmesh_error_msg("i%2 must be either 0 or 1!");
265  }
266 
267  // dummy
268  return libMesh::RealVectorValue();
269 }
270 
271 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
272 template <>
275  const Order order,
276  const unsigned int i,
277  const unsigned int j,
278  const Point & p)
279 {
280  Real value = FE<2, MONOMIAL>::shape_second_deriv(type, order, i / 2, j, p);
281 
282  switch (i % 2)
283  {
284  case 0:
286 
287  case 1:
288  return libMesh::RealVectorValue(Real(0), value);
289 
290  default:
291  libmesh_error_msg("i%2 must be either 0 or 1!");
292  }
293 
294  // dummy
295  return libMesh::RealVectorValue();
296 }
297 
298 #endif
299 
300 // 3-D
301 template <>
304  const Order order,
305  const unsigned int i,
306  const Point & p)
307 {
308  Real value = FE<3, MONOMIAL>::shape(type, order, i / 3, p);
309 
310  switch (i % 3)
311  {
312  case 0:
314 
315  case 1:
316  return libMesh::RealVectorValue(Real(0), value);
317 
318  case 2:
319  return libMesh::RealVectorValue(Real(0), Real(0), value);
320 
321  default:
322  libmesh_error_msg("i%3 must be 0, 1, or 2!");
323  }
324 
325  // dummy
326  return libMesh::RealVectorValue();
327 }
328 template <>
331  const Order order,
332  const unsigned int i,
333  const unsigned int j,
334  const Point & p)
335 {
336  Real value = FE<3, MONOMIAL>::shape_deriv(type, order, i / 3, j, p);
337 
338  switch (i % 3)
339  {
340  case 0:
342 
343  case 1:
344  return libMesh::RealVectorValue(Real(0), value);
345 
346  case 2:
347  return libMesh::RealVectorValue(Real(0), Real(0), value);
348 
349  default:
350  libmesh_error_msg("i%3 must be 0, 1, or 2!");
351  }
352 
353  // dummy
354  return libMesh::RealVectorValue();
355 }
356 
357 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
358 
359 template <>
362  const Order order,
363  const unsigned int i,
364  const unsigned int j,
365  const Point & p)
366 {
367  Real value = FE<3, MONOMIAL>::shape_second_deriv(type, order, i / 3, j, p);
368 
369  switch (i % 3)
370  {
371  case 0:
373 
374  case 1:
375  return libMesh::RealVectorValue(Real(0), value);
376 
377  case 2:
378  return libMesh::RealVectorValue(Real(0), Real(0), value);
379 
380  default:
381  libmesh_error_msg("i%3 must be 0, 1, or 2!");
382  }
383 
384  // dummy
385  return libMesh::RealVectorValue();
386 }
387 
388 #endif
389 
390 // 0-D
391 template <>
394  const Order order,
395  const unsigned int i,
396  const Point & p,
397  const bool add_p_level)
398 {
399  Real value =
400  FE<0, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i, p);
402 }
403 template <>
406  const Order order,
407  const unsigned int i,
408  const unsigned int j,
409  const Point & p,
410  const bool add_p_level)
411 {
413  elem->type(), order + add_p_level*elem->p_level(), i, j, p);
415 }
416 
417 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
418 
419 template <>
422  const Order order,
423  const unsigned int i,
424  const unsigned int j,
425  const Point & p,
426  const bool add_p_level)
427 {
429  elem->type(), order + add_p_level*elem->p_level(), i, j, p);
431 }
432 
433 #endif
434 
435 // 1-D
436 template <>
439  const Order order,
440  const unsigned int i,
441  const Point & p,
442  const bool add_p_level)
443 {
444  Real value =
445  FE<1, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i, p);
447 }
448 template <>
451  const Order order,
452  const unsigned int i,
453  const unsigned int j,
454  const Point & p,
455  const bool add_p_level)
456 {
458  elem->type(), order + add_p_level*elem->p_level(), i, j, p);
460 }
461 
462 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
463 template <>
466  const Order order,
467  const unsigned int i,
468  const unsigned int j,
469  const Point & p,
470  const bool add_p_level)
471 {
473  elem->type(), order + add_p_level*elem->p_level(), i, j, p);
475 }
476 
477 #endif
478 
479 // 2-D
480 template <>
483  const Order order,
484  const unsigned int i,
485  const Point & p,
486  const bool add_p_level)
487 {
488  Real value =
489  FE<2, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i / 2, p);
490 
491  switch (i % 2)
492  {
493  case 0:
495 
496  case 1:
497  return libMesh::RealVectorValue(Real(0), value);
498 
499  default:
500  libmesh_error_msg("i%2 must be either 0 or 1!");
501  }
502 
503  // dummy
504  return libMesh::RealVectorValue();
505 }
506 template <>
509  const Order order,
510  const unsigned int i,
511  const unsigned int j,
512  const Point & p,
513  const bool add_p_level)
514 {
516  elem->type(), order + add_p_level*elem->p_level(), i / 2, j, p);
517 
518  switch (i % 2)
519  {
520  case 0:
522 
523  case 1:
524  return libMesh::RealVectorValue(Real(0), value);
525 
526  default:
527  libmesh_error_msg("i%2 must be either 0 or 1!");
528  }
529 
530  // dummy
531  return libMesh::RealVectorValue();
532 }
533 
534 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
535 template <>
538  const Order order,
539  const unsigned int i,
540  const unsigned int j,
541  const Point & p,
542  const bool add_p_level)
543 {
545  elem->type(), order + add_p_level*elem->p_level(), i / 2, j, p);
546 
547  switch (i % 2)
548  {
549  case 0:
551 
552  case 1:
553  return libMesh::RealVectorValue(Real(0), value);
554 
555  default:
556  libmesh_error_msg("i%2 must be either 0 or 1!");
557  }
558 
559  // dummy
560  return libMesh::RealVectorValue();
561 }
562 
563 #endif
564 
565 // 3-D
566 template <>
569  const Order order,
570  const unsigned int i,
571  const Point & p,
572  const bool add_p_level)
573 {
574  Real value =
575  FE<3, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i / 3, p);
576 
577  switch (i % 3)
578  {
579  case 0:
581 
582  case 1:
583  return libMesh::RealVectorValue(Real(0), value);
584 
585  case 2:
586  return libMesh::RealVectorValue(Real(0), Real(0), value);
587 
588  default:
589  libmesh_error_msg("i%3 must be 0, 1, or 2!");
590  }
591 
592  // dummy
593  return libMesh::RealVectorValue();
594 }
595 template <>
598  const Order order,
599  const unsigned int i,
600  const unsigned int j,
601  const Point & p,
602  const bool add_p_level)
603 {
605  elem->type(), order + add_p_level*elem->p_level(), i / 3, j, p);
606 
607  switch (i % 3)
608  {
609  case 0:
611 
612  case 1:
613  return libMesh::RealVectorValue(Real(0), value);
614 
615  case 2:
616  return libMesh::RealVectorValue(Real(0), Real(0), value);
617 
618  default:
619  libmesh_error_msg("i%3 must be 0, 1, or 2!");
620  }
621 
622  // dummy
623  return libMesh::RealVectorValue();
624 }
625 
626 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
627 
628 template <>
631  const Order order,
632  const unsigned int i,
633  const unsigned int j,
634  const Point & p,
635  const bool add_p_level)
636 {
638  elem->type(), order + add_p_level*elem->p_level(), i / 3, j, p);
639 
640  switch (i % 3)
641  {
642  case 0:
644 
645  case 1:
646  return libMesh::RealVectorValue(Real(0), value);
647 
648  case 2:
649  return libMesh::RealVectorValue(Real(0), Real(0), value);
650 
651  default:
652  libmesh_error_msg("i%3 must be 0, 1, or 2!");
653  }
654 
655  // dummy
656  return libMesh::RealVectorValue();
657 }
658 
659 #endif
660 
661 // Do full-specialization for every dimension, instead
662 // of explicit instantiation at the end of this function.
663 // This could be macro-ified.
664 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0, MONOMIAL>::n_dofs(t, o); }
665 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1, MONOMIAL>::n_dofs(t, o); }
666 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs(t, o); }
667 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs(t, o); }
668 
669 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0, MONOMIAL>::n_dofs(e, o); }
670 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1, MONOMIAL>::n_dofs(e, o); }
671 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs(e, o); }
672 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs(e, o); }
673 
674 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
675 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
676 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
677 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
678 
679 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
680 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
681 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
682 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
683 
684 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0, MONOMIAL>::n_dofs_per_elem(t, o); }
685 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1, MONOMIAL>::n_dofs_per_elem(t, o); }
686 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs_per_elem(t, o); }
687 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs_per_elem(t, o); }
688 
689 template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0, MONOMIAL>::n_dofs_per_elem(e, o); }
690 template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1, MONOMIAL>::n_dofs_per_elem(e, o); }
691 template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs_per_elem(e, o); }
692 template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs_per_elem(e, o); }
693 
694 // Monomial FEMs are always C^0 continuous
695 template <>
698 {
699  return DISCONTINUOUS;
700 }
701 template <>
704 {
705  return DISCONTINUOUS;
706 }
707 template <>
710 {
711  return DISCONTINUOUS;
712 }
713 template <>
716 {
717  return DISCONTINUOUS;
718 }
719 
720 // Monomial FEMs are not hierarchic
721 template <>
722 bool
724 {
725  return true;
726 }
727 template <>
728 bool
730 {
731  return true;
732 }
733 template <>
734 bool
736 {
737  return true;
738 }
739 template <>
740 bool
742 {
743  return true;
744 }
745 
746 // Monomial FEM shapes do not need reinit (is this always true?)
747 template <>
748 bool
750 {
751  return false;
752 }
753 template <>
754 bool
756 {
757  return false;
758 }
759 template <>
760 bool
762 {
763  return false;
764 }
765 template <>
766 bool
768 {
769  return false;
770 }
771 
772 // we only need instantiations of this function for
773 // Dim==2 and 3. There are no constraints for discontinuous elements, so
774 // we do nothing.
775 #ifdef LIBMESH_ENABLE_AMR
776 template <>
777 void
779  DofMap &,
780  const unsigned int,
781  const Elem *)
782 {
783 }
784 
785 template <>
786 void
788  DofMap &,
789  const unsigned int,
790  const Elem *)
791 {
792 }
793 #endif // LIBMESH_ENABLE_AMR
794 
795 } // namespace libMesh
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
static unsigned int n_dofs(const ElemType t, const Order o)
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
LIBMESH_FE_NODAL_SOLN_DIM(LIBMESH_FE_NODAL_SOLN_DIM(HIERARCHIC_VEC,(FE< 0, HIERARCHIC >::nodal_soln), 0)
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
unsigned int dim
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
Definition: elem.h:3108
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
A specific instantiation of the FEBase class.
Definition: fe.h:127
LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:272
virtual unsigned int n_nodes() const =0
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:760
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:400
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln from the element soln.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
static const bool value
Definition: xdr_io.C:54
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
The constraint matrix storage format.
Definition: dof_map.h:108
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)