Intrepid2
Intrepid2_DerivedBasis_HCURL_HEX.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
20 #ifndef Intrepid2_DerivedBasis_HCURL_HEX_h
21 #define Intrepid2_DerivedBasis_HCURL_HEX_h
22 
23 #include <Kokkos_DynRankView.hpp>
24 
26 #include "Intrepid2_Sacado.hpp"
27 
30 
31 namespace Intrepid2
32 {
33  template<class HGRAD_LINE, class HVOL_LINE>
35  : public Basis_TensorBasis3<typename HGRAD_LINE::BasisBase>
36  {
37  public:
38  using OutputViewType = typename HGRAD_LINE::OutputViewType;
39  using PointViewType = typename HGRAD_LINE::PointViewType ;
40  using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
41 
42  using LineGradBasis = HGRAD_LINE;
43  using LineVolBasis = HVOL_LINE;
44 
46  public:
53  Basis_Derived_HCURL_Family1_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType = POINTTYPE_DEFAULT)
54  :
55  TensorBasis3(Teuchos::rcp(new LineVolBasis (polyOrder_x-1,pointType)),
56  Teuchos::rcp(new LineGradBasis(polyOrder_y,pointType)),
57  Teuchos::rcp(new LineGradBasis(polyOrder_z,pointType)),
58  true) // true: use shards CellTopology and tags
59  {
60  this->functionSpace_ = FUNCTION_SPACE_HCURL;
61  }
62 
65  virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
66  {
67  const EOperator VALUE = Intrepid2::OPERATOR_VALUE;
68  const EOperator GRAD = Intrepid2::OPERATOR_GRAD;
69 
70  if (operatorType == Intrepid2::OPERATOR_VALUE)
71  {
72  // family 1 goes in x component
73  std::vector< std::vector<EOperator> > ops(3);
74  ops[0] = std::vector<EOperator>{VALUE,VALUE,VALUE};
75  ops[1] = std::vector<EOperator>{};
76  ops[2] = std::vector<EOperator>{};
77  std::vector<double> weights {1.0, 0.0, 0.0};
78  return OperatorTensorDecomposition(ops, weights);
79  }
80  else if (operatorType == Intrepid2::OPERATOR_CURL)
81  {
82  // Family 1:
83  // x component is zero
84  // y component is d/dz: (VALUE,VALUE,GRAD), weight = 1.0
85  // z component is -d/dy: (VALUE,GRAD,VALUE), weight = -1.0
86 
87  std::vector< std::vector<EOperator> > ops(3);
88  ops[0] = std::vector<EOperator>{};
89  ops[1] = std::vector<EOperator>{VALUE,VALUE,GRAD};
90  ops[2] = std::vector<EOperator>{VALUE,GRAD,VALUE};
91 
92  std::vector<double> weights {0.0, 1.0, -1.0};
93  return OperatorTensorDecomposition(ops,weights);
94  }
95  else
96  {
97  INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported operator type");
98  }
99  }
100 
102 
111  virtual void getValues(OutputViewType outputValues, const EOperator operatorType,
112  const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3,
113  bool tensorPoints) const override
114  {
115  Intrepid2::EOperator op1, op2, op3;
116  if (operatorType == Intrepid2::OPERATOR_VALUE)
117  {
118  op1 = Intrepid2::OPERATOR_VALUE;
119  op2 = Intrepid2::OPERATOR_VALUE;
120  op3 = Intrepid2::OPERATOR_VALUE;
121 
122  // family 1 goes in the x component; 0 in the y and z components
123  auto outputValuesComponent1 = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),0);
124  auto outputValuesComponent23 = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),std::make_pair(1,3));
125 
126  this->TensorBasis3::getValues(outputValuesComponent1,
127  inputPoints1, op1,
128  inputPoints2, op2,
129  inputPoints3, op3, tensorPoints);
130  // place 0 in the y and z components
131  Kokkos::deep_copy(outputValuesComponent23,0.0);
132  }
133  else if (operatorType == Intrepid2::OPERATOR_CURL)
134  {
135  // family 1 is nonzero in the x component, so the curl is d/dz placed in the y component, and -d/dy placed in the z component.
136  auto outputValuesComponent_x = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),0);
137  auto outputValuesComponent_y = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),1);
138  auto outputValuesComponent_z = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),2);
139 
140  // x component is zero
141  Kokkos::deep_copy(outputValuesComponent_x, 0.0);
142 
143  // y component is d/dz
144  op1 = Intrepid2::OPERATOR_VALUE;
145  op2 = Intrepid2::OPERATOR_VALUE;
146  op3 = Intrepid2::OPERATOR_GRAD; // d/dz
147 
148  double weight = 1.0; // the plus sign in front of d/dz
149  this->TensorBasis3::getValues(outputValuesComponent_y,
150  inputPoints1, op1,
151  inputPoints2, op2,
152  inputPoints3, op3, tensorPoints, weight);
153 
154  // z component is -d/dy
155  op1 = Intrepid2::OPERATOR_VALUE;
156  op2 = Intrepid2::OPERATOR_GRAD; // d/dy
157  op3 = Intrepid2::OPERATOR_VALUE;
158  weight = -1.0; // the -1 weight on d/dy
159  this->TensorBasis3::getValues(outputValuesComponent_z,
160  inputPoints1, op1,
161  inputPoints2, op2,
162  inputPoints3, op3, tensorPoints, weight);
163  }
164  else
165  {
166  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"operator not yet supported");
167  }
168  }
169 
181  virtual void getDofCoeffs( ScalarViewType dofCoeffs ) const override {
182  auto dofCoeffs1 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),0);
183  auto dofCoeffs23 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),std::make_pair(1,3));
184  this->TensorBasis3::getDofCoeffs(dofCoeffs1);
185  Kokkos::deep_copy(dofCoeffs23,0.0);
186  }
187  };
188 
189  template<class HGRAD_LINE, class HVOL_LINE>
191  : public Basis_TensorBasis3<typename HGRAD_LINE::BasisBase>
192  {
193  public:
194  using OutputViewType = typename HGRAD_LINE::OutputViewType;
195  using PointViewType = typename HGRAD_LINE::PointViewType ;
196  using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
197 
198  using LineGradBasis = HGRAD_LINE;
199  using LineVolBasis = HVOL_LINE;
200 
202  public:
209  Basis_Derived_HCURL_Family2_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType = POINTTYPE_DEFAULT)
210  :
211  TensorBasis3(Teuchos::rcp( new LineGradBasis(polyOrder_x,pointType)),
212  Teuchos::rcp( new LineVolBasis (polyOrder_y-1,pointType)),
213  Teuchos::rcp( new LineGradBasis(polyOrder_z,pointType)),
214  true) // true: use shards CellTopology and tags
215  {
216  this->functionSpace_ = FUNCTION_SPACE_HCURL;
217  }
218 
221  virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
222  {
223  const EOperator VALUE = Intrepid2::OPERATOR_VALUE;
224  const EOperator GRAD = Intrepid2::OPERATOR_GRAD;
225  const EOperator CURL = Intrepid2::OPERATOR_CURL;
226  if (operatorType == VALUE)
227  {
228  // family 2 goes in y component
229  std::vector< std::vector<EOperator> > ops(3);
230  ops[0] = std::vector<EOperator>{};
231  ops[1] = std::vector<EOperator>{VALUE,VALUE,VALUE};
232  ops[2] = std::vector<EOperator>{};
233  std::vector<double> weights {0.0, 1.0, 0.0};
234  return OperatorTensorDecomposition(ops, weights);
235  }
236  else if (operatorType == CURL)
237  {
238  // family 2 is nonzero in the y component, so the curl is -d/dz placed in the x component, and d/dx placed in the z component.
239  // x component is -d/dz: (VALUE,VALUE,GRAD), weight = -1.0
240  // y component is zero
241  // z component is d/dx: (GRAD,VALUE,VALUE), weight = 1.0
242  std::vector< std::vector<EOperator> > ops(3);
243  ops[0] = std::vector<EOperator>{VALUE,VALUE,GRAD};
244  ops[1] = std::vector<EOperator>{};
245  ops[2] = std::vector<EOperator>{GRAD,VALUE,VALUE};
246 
247  std::vector<double> weights {-1.0, 0.0, 1.0};
248  return OperatorTensorDecomposition(ops,weights);
249  }
250  else
251  {
252  INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported operator type");
253  }
254  }
255 
257 
266  virtual void getValues(OutputViewType outputValues, const EOperator operatorType,
267  const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3,
268  bool tensorPoints) const override
269  {
270  Intrepid2::EOperator op1, op2, op3;
271  if (operatorType == Intrepid2::OPERATOR_VALUE)
272  {
273  op1 = Intrepid2::OPERATOR_VALUE;
274  op2 = Intrepid2::OPERATOR_VALUE;
275  op3 = Intrepid2::OPERATOR_VALUE;
276 
277  // family 2 goes in the y component; 0 in the x and z components
278  auto outputValuesComponent_x = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),0);
279  auto outputValuesComponent_y = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),1);
280  auto outputValuesComponent_z = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),2);
281 
282  // place 0 in the x component
283  Kokkos::deep_copy(outputValuesComponent_x,0.0);
284  // evaluate y component
285  this->TensorBasis3::getValues(outputValuesComponent_y,
286  inputPoints1, op1,
287  inputPoints2, op2,
288  inputPoints3, op3, tensorPoints);
289  // place 0 in the z component
290  Kokkos::deep_copy(outputValuesComponent_z,0.0);
291  }
292  else if (operatorType == Intrepid2::OPERATOR_CURL)
293  {
294  // family 2 is nonzero in the y component, so the curl is -d/dz placed in the x component, and d/dx placed in the z component.
295  auto outputValuesComponent_x = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),0);
296  auto outputValuesComponent_y = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),1);
297  auto outputValuesComponent_z = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),2);
298 
299  // x component is -d/dz
300  op1 = Intrepid2::OPERATOR_VALUE;
301  op2 = Intrepid2::OPERATOR_VALUE;
302  op3 = Intrepid2::OPERATOR_GRAD; // d/dz
303 
304  double weight = -1.0; // the minus sign in front of d/dz
305  this->TensorBasis3::getValues(outputValuesComponent_x,
306  inputPoints1, op1,
307  inputPoints2, op2,
308  inputPoints3, op3, tensorPoints, weight);
309 
310  // y component is zero
311  Kokkos::deep_copy(outputValuesComponent_y, 0.0);
312 
313  // z component is d/dx
314  op1 = Intrepid2::OPERATOR_GRAD; // d/dx
315  op2 = Intrepid2::OPERATOR_VALUE;
316  op3 = Intrepid2::OPERATOR_VALUE;
317  weight = 1.0; // the weight on d/dx
318  this->TensorBasis3::getValues(outputValuesComponent_z,
319  inputPoints1, op1,
320  inputPoints2, op2,
321  inputPoints3, op3, tensorPoints, weight);
322  }
323  else
324  {
325  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"operator not yet supported");
326  }
327  }
328 
340  virtual void getDofCoeffs( ScalarViewType dofCoeffs ) const override {
341  auto dofCoeffs1 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),0);
342  auto dofCoeffs2 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),1);
343  auto dofCoeffs3 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),2);
344  Kokkos::deep_copy(dofCoeffs1,0.0);
345  this->TensorBasis3::getDofCoeffs(dofCoeffs2);
346  Kokkos::deep_copy(dofCoeffs3,0.0);
347  }
348 
349 
350  };
351 
352  template<class HGRAD_LINE, class HVOL_LINE>
354  : public Basis_TensorBasis3<typename HGRAD_LINE::BasisBase>
355  {
356  using OutputViewType = typename HGRAD_LINE::OutputViewType;
357  using PointViewType = typename HGRAD_LINE::PointViewType ;
358  using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
359 
360  using LineGradBasis = HGRAD_LINE;
361  using LineVolBasis = HVOL_LINE;
362 
364  public:
371  Basis_Derived_HCURL_Family3_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType = POINTTYPE_DEFAULT)
372  :
373  TensorBasis3(Teuchos::rcp(new LineGradBasis(polyOrder_x,pointType)),
374  Teuchos::rcp(new LineGradBasis(polyOrder_y,pointType)),
375  Teuchos::rcp(new LineVolBasis (polyOrder_z-1,pointType)),
376  true) // true: use shards CellTopology and tags
377  {
378  this->functionSpace_ = FUNCTION_SPACE_HCURL;
379  }
380 
383  virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
384  {
385  const EOperator VALUE = Intrepid2::OPERATOR_VALUE;
386  const EOperator GRAD = Intrepid2::OPERATOR_GRAD;
387  if (operatorType == Intrepid2::OPERATOR_VALUE)
388  {
389  // family 3 goes in z component
390  std::vector< std::vector<EOperator> > ops(3);
391  ops[0] = std::vector<EOperator>{};
392  ops[1] = std::vector<EOperator>{};
393  ops[2] = std::vector<EOperator>{VALUE,VALUE,VALUE};
394  std::vector<double> weights {0.0, 0.0, 1.0};
395  return OperatorTensorDecomposition(ops, weights);
396  }
397  else if (operatorType == Intrepid2::OPERATOR_CURL)
398  {
399  // family 3 is nonzero in the z component, so the curl is d/dy placed in the x component, and -d/dx placed in the y component.
400  // x component is d/dy: (VALUE,GRAD,VALUE), weight = 1.0
401  // y component is d/dx: (GRAD,VALUE,VALUE), weight = -1.0
402  // z component is zero
403  std::vector< std::vector<EOperator> > ops(3);
404  ops[0] = std::vector<EOperator>{VALUE,GRAD,VALUE};
405  ops[1] = std::vector<EOperator>{GRAD,VALUE,VALUE};
406  ops[2] = std::vector<EOperator>{};
407 
408  std::vector<double> weights {1.0, -1.0, 0.0};
409  return OperatorTensorDecomposition(ops,weights);
410  }
411  else
412  {
413  INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported operator type");
414  }
415  }
416 
418 
427  virtual void getValues(OutputViewType outputValues, const EOperator operatorType,
428  const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3,
429  bool tensorPoints) const override
430  {
431  Intrepid2::EOperator op1, op2, op3;
432  if (operatorType == Intrepid2::OPERATOR_VALUE)
433  {
434  op1 = Intrepid2::OPERATOR_VALUE;
435  op2 = Intrepid2::OPERATOR_VALUE;
436  op3 = Intrepid2::OPERATOR_VALUE;
437 
438  // family 3 goes in the z component; 0 in the x and y components
439  auto outputValuesComponent_xy = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),std::make_pair(0,2));
440  auto outputValuesComponent_z = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),2);
441 
442  // place 0 in the x and y components
443  Kokkos::deep_copy(outputValuesComponent_xy,0.0);
444  // evaluate z component
445  this->TensorBasis3::getValues(outputValuesComponent_z,
446  inputPoints1, op1,
447  inputPoints2, op2,
448  inputPoints3, op3, tensorPoints);
449  }
450  else if (operatorType == Intrepid2::OPERATOR_CURL)
451  {
452  // family 3 is nonzero in the z component, so the curl is d/dy placed in the x component, and -d/dx placed in the y component.
453  auto outputValuesComponent_x = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),0);
454  auto outputValuesComponent_y = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),1);
455  auto outputValuesComponent_z = Kokkos::subview(outputValues,Kokkos::ALL(),Kokkos::ALL(),2);
456 
457  // x component is d/dy
458  op1 = Intrepid2::OPERATOR_VALUE;
459  op2 = Intrepid2::OPERATOR_GRAD; // d/dy
460  op3 = Intrepid2::OPERATOR_VALUE;
461 
462  double weight = 1.0; // the sign in front of d/dy
463  this->TensorBasis3::getValues(outputValuesComponent_x,
464  inputPoints1, op1,
465  inputPoints2, op2,
466  inputPoints3, op3, tensorPoints, weight);
467  // y component is -d/dx
468  op1 = Intrepid2::OPERATOR_GRAD; // d/dx
469  op2 = Intrepid2::OPERATOR_VALUE;
470  op3 = Intrepid2::OPERATOR_VALUE;
471  weight = -1.0; // the weight on d/dx
472  this->TensorBasis3::getValues(outputValuesComponent_y,
473  inputPoints1, op1,
474  inputPoints2, op2,
475  inputPoints3, op3, tensorPoints, weight);
476 
477  // z component is zero
478  Kokkos::deep_copy(outputValuesComponent_z, 0.0);
479  }
480  else
481  {
482  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"operator not yet supported");
483  }
484  }
485 
497  virtual void getDofCoeffs( ScalarViewType dofCoeffs ) const override {
498  auto dofCoeffs12 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),std::make_pair(0,2));
499  auto dofCoeffs3 = Kokkos::subview(dofCoeffs,Kokkos::ALL(),2);
500  Kokkos::deep_copy(dofCoeffs12,0.0);
501  this->TensorBasis3::getDofCoeffs(dofCoeffs3);
502  }
503  };
504 
505  template<class HGRAD_LINE, class HVOL_LINE>
507  : public Basis_DirectSumBasis <typename HGRAD_LINE::BasisBase>
508  {
512  public:
519  Basis_Derived_HCURL_Family1_Family2_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType)
520  :
521  DirectSumBasis(Teuchos::rcp(new Family1(polyOrder_x, polyOrder_y, polyOrder_z, pointType)),
522  Teuchos::rcp(new Family2(polyOrder_x, polyOrder_y, polyOrder_z, pointType))) {}
523  };
524 
525  template<class HGRAD_LINE, class HVOL_LINE>
527  : public Basis_DirectSumBasis <typename HGRAD_LINE::BasisBase>
528  {
532 
533  std::string name_;
534  ordinal_type order_x_;
535  ordinal_type order_y_;
536  ordinal_type order_z_;
537  EPointType pointType_;
538 
539  public:
540 
541  using ExecutionSpace = typename HGRAD_LINE::ExecutionSpace;
542  using OutputValueType = typename HGRAD_LINE::OutputValueType;
543  using PointValueType = typename HGRAD_LINE::PointValueType;
544 
545  using BasisBase = typename HGRAD_LINE::BasisBase;
546 
553  Basis_Derived_HCURL_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
554  :
555  DirectSumBasis(Teuchos::rcp(new Family12(polyOrder_x, polyOrder_y, polyOrder_z, pointType)),
556  Teuchos::rcp(new Family3 (polyOrder_x, polyOrder_y, polyOrder_z, pointType))) {
557  this->functionSpace_ = FUNCTION_SPACE_HCURL;
558 
559  std::ostringstream basisName;
560  basisName << "HCURL_HEX (" << this->DirectSumBasis::getName() << ")";
561  name_ = basisName.str();
562 
563  order_x_ = polyOrder_x;
564  order_y_ = polyOrder_y;
565  order_z_ = polyOrder_z;
566  pointType_ = pointType;
567  }
568 
572  Basis_Derived_HCURL_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) : Basis_Derived_HCURL_HEX(polyOrder, polyOrder, polyOrder, pointType) {}
573 
576  virtual bool requireOrientation() const override {
577  return true;
578  }
579 
584  virtual
585  const char*
586  getName() const override {
587  return name_.c_str();
588  }
589 
600  Teuchos::RCP<BasisBase>
601  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
602 
603  using LineBasis = HVOL_LINE;
605 
606  if(subCellDim == 1) {
607  switch(subCellOrd) {
608  case 0:
609  case 2:
610  case 4:
611  case 6:
612  return Teuchos::rcp( new LineBasis(order_x_-1, pointType_) );
613  case 1:
614  case 3:
615  case 5:
616  case 7:
617  return Teuchos::rcp( new LineBasis(order_y_-1, pointType_) );
618  case 8:
619  case 9:
620  case 10:
621  case 11:
622  return Teuchos::rcp( new LineBasis(order_z_-1, pointType_) );
623  }
624  } else if(subCellDim == 2) {
625  switch(subCellOrd) {
626  case 0:
627  return Teuchos::rcp( new QuadBasis(order_x_, order_z_, pointType_) );
628  case 1:
629  return Teuchos::rcp( new QuadBasis(order_y_,order_z_, pointType_) );
630  case 2:
631  return Teuchos::rcp( new QuadBasis(order_x_, order_z_, pointType_) );
632  case 3:
633  return Teuchos::rcp( new QuadBasis(order_z_, order_y_, pointType_) );
634  case 4:
635  return Teuchos::rcp( new QuadBasis(order_y_, order_x_, pointType_) );
636  case 5:
637  return Teuchos::rcp( new QuadBasis(order_x_, order_y_, pointType_) );
638  }
639  }
640 
641  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
642  }
647  virtual HostBasisPtr<OutputValueType, PointValueType>
648  getHostBasis() const override {
650 
651  auto hostBasis = Teuchos::rcp(new HostBasis(order_x_, order_y_, order_z_, pointType_));
652 
653  return hostBasis;
654  }
655  };
656 } // end namespace Intrepid2
657 
658 #endif /* Intrepid2_DerivedBasis_HCURL_HEX_h */
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Fills in coefficients of degrees of freedom for Lagrangian basis on the reference cell...
Teuchos::RCP< BasisBase > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual bool requireOrientation() const override
True if orientation is required.
Basis_Derived_HCURL_Family1_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis_Derived_HCURL_Family1_Family2_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType)
Constructor.
virtual const char * getName() const override
Returns basis name.
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis...
Implementation of bases that are tensor products of two or three component bases. ...
Free functions, callable from device code, that implement various polynomials useful in basis definit...
For a multi-component tensor basis, specifies the operators to be applied to the components to produc...
virtual void getValues(OutputViewType outputValues, const EOperator operatorType, const PointViewType inputPoints12, const PointViewType inputPoints3, bool tensorPoints) const override
Evaluation of a tensor FEM basis on a reference cell.
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis...
Basis_Derived_HCURL_Family2_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Implementation of a basis that is the direct sum of two other bases.
A basis that is the direct sum of two other bases.
virtual const char * getName() const override
Returns basis name.
Basis_Derived_HCURL_Family3_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis_Derived_HCURL_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual void getValues(OutputViewType outputValues, const EOperator operatorType, const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3, bool tensorPoints) const override
multi-component getValues() method (required/called by TensorBasis3)
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Fills in coefficients of degrees of freedom for Lagrangian basis on the reference cell...
virtual void getValues(OutputViewType outputValues, const EOperator operatorType, const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3, bool tensorPoints) const override
multi-component getValues() method (required/called by TensorBasis3)
Header file to include all Sacado headers that are required if using Intrepid2 with Sacado types...
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis...
virtual void getDofCoeffs(typename BasisBase::ScalarViewType dofCoeffs) const override
Fills in coefficients of degrees of freedom on the reference cell.
virtual void getValues(OutputViewType outputValues, const EOperator operatorType, const PointViewType inputPoints1, const PointViewType inputPoints2, const PointViewType inputPoints3, bool tensorPoints) const override
multi-component getValues() method (required/called by TensorBasis3)
Basis_Derived_HCURL_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual HostBasisPtr< OutputValueType, PointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Fills in coefficients of degrees of freedom for Lagrangian basis on the reference cell...