Intrepid2
Intrepid2_HDIV_TET_In_FEMDef.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
49#ifndef __INTREPID2_HDIV_TET_IN_FEM_DEF_HPP__
50#define __INTREPID2_HDIV_TET_IN_FEM_DEF_HPP__
51
54
55namespace Intrepid2 {
56
57// -------------------------------------------------------------------------------------
58
59namespace Impl {
60
61template<EOperator opType>
62template<typename OutputViewType,
63typename inputViewType,
64typename workViewType,
65typename vinvViewType>
66KOKKOS_INLINE_FUNCTION
67void
68Basis_HDIV_TET_In_FEM::Serial<opType>::
69getValues( OutputViewType output,
70 const inputViewType input,
71 workViewType work,
72 const vinvViewType coeffs ) {
73
74 constexpr ordinal_type spaceDim = 3;
75 const ordinal_type
76 cardPn = coeffs.extent(0)/spaceDim,
77 card = coeffs.extent(1),
78 npts = input.extent(0);
79
80 // compute order
81 ordinal_type order = 0;
82 for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
83 if (card == CardinalityHDivTet(p)) {
84 order = p;
85 break;
86 }
87 }
88
89 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
90 auto vcprop = Kokkos::common_view_alloc_prop(work);
91 auto ptr = work.data();
92
93 switch (opType) {
94 case OPERATOR_VALUE: {
95 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
96 workViewType dummyView;
97
98 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
99 Serial<opType>::getValues(phis, input, dummyView, order);
100
101 for (ordinal_type i=0;i<card;++i)
102 for (ordinal_type j=0;j<npts;++j)
103 for (ordinal_type d=0;d<spaceDim;++d) {
104 output.access(i,j,d) = 0.0;
105 for (ordinal_type k=0;k<cardPn;++k)
106 output.access(i,j,d) += coeffs(k+d*cardPn,i) * phis.access(k,j);
107 }
108 break;
109 }
110 case OPERATOR_DIV: {
111 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
112 ptr += card*npts*spaceDim*get_dimension_scalar(work);
113 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
114
115 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
116 Serial<OPERATOR_GRAD>::getValues(phis, input, workView, order);
117
118 for (ordinal_type i=0;i<card;++i)
119 for (ordinal_type j=0;j<npts;++j) {
120 output.access(i,j) = 0.0;
121 for (ordinal_type k=0; k<cardPn; ++k)
122 for (ordinal_type d=0; d<spaceDim; ++d)
123 output.access(i,j) += coeffs(k+d*cardPn,i)*phis.access(k,j,d);
124 }
125 break;
126 }
127 default: {
128 INTREPID2_TEST_FOR_ABORT( true,
129 ">>> ERROR (Basis_HDIV_TET_In_FEM): Operator type not implemented");
130 }
131 }
132}
133
134template<typename DT, ordinal_type numPtsPerEval,
135typename outputValueValueType, class ...outputValueProperties,
136typename inputPointValueType, class ...inputPointProperties,
137typename vinvValueType, class ...vinvProperties>
138void
139Basis_HDIV_TET_In_FEM::
140getValues( /* */ Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
141 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
142 const Kokkos::DynRankView<vinvValueType, vinvProperties...> coeffs,
143 const EOperator operatorType) {
144 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
145 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
146 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
147 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
148
149 // loopSize corresponds to cardinality
150 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
151 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
152 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
153 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
154
155 typedef typename inputPointViewType::value_type inputPointType;
156
157 const ordinal_type cardinality = outputValues.extent(0);
158 const ordinal_type spaceDim = 3;
159
160 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
161 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
162
163 switch (operatorType) {
164 case OPERATOR_VALUE: {
165 workViewType work(Kokkos::view_alloc("Basis_HDIV_TET_In_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
166 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
167 OPERATOR_VALUE,numPtsPerEval> FunctorType;
168 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
169 break;
170 }
171 case OPERATOR_DIV: {
172 workViewType work(Kokkos::view_alloc("Basis_HDIV_TET_In_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
173 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
174 OPERATOR_DIV,numPtsPerEval> FunctorType;
175 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
176 break;
177 }
178 default: {
179 INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
180 ">>> ERROR (Basis_HDIV_TET_In_FEM): Operator type not implemented" );
181 }
182 }
183}
184}
185
186// -------------------------------------------------------------------------------------
187template<typename DT, typename OT, typename PT>
189Basis_HDIV_TET_In_FEM( const ordinal_type order,
190 const EPointType pointType ) {
191
192 constexpr ordinal_type spaceDim = 3;
193 this->basisCardinality_ = CardinalityHDivTet(order);
194 this->basisDegree_ = order; // small n
195 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
196 this->basisType_ = BASIS_FEM_LAGRANGIAN;
197 this->basisCoordinates_ = COORDINATES_CARTESIAN;
198 this->functionSpace_ = FUNCTION_SPACE_HDIV;
199 pointType_ = pointType;
200
201 const ordinal_type card = this->basisCardinality_;
202
203 const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order); // dim of (P_{n}) -- smaller space
204 const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1); // dim of (P_{n-1}) -- smaller space
205 const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2); // dim of (P_{n-2}) -- smaller space
206 const ordinal_type cardVecPn = spaceDim*cardPn; // dim of (P_{n})^3 -- larger space
207 const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^3 -- smaller space
208 const ordinal_type dim_PkH = cardPnm1 - cardPnm2;
209
210
211 // Basis-dependent initializations
212 constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
213 constexpr ordinal_type maxCard = CardinalityHDivTet(Parameters::MaxOrder);
214 ordinal_type tags[maxCard][tagSize];
215
216 // points are computed in the host and will be copied
217 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
218 dofCoords("Hdiv::Tet::In::dofCoords", card, spaceDim);
219
220 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
221 dofCoeffs("Hdiv::Tet::In::dofCoeffs", card, spaceDim);
222
223 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
224 coeffs("Hdiv::Tet::In::coeffs", cardVecPn, card);
225
226 // first, need to project the basis for RT space onto the
227 // orthogonal basis of degree n
228 // get coefficients of PkHx
229
230 const ordinal_type lwork = card*card;
231 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
232 V1("Hdiv::Tet::In::V1", cardVecPn, card);
233
234 // basis for the space is
235 // { (phi_i,0,0) }_{i=0}^{cardPnm1-1} ,
236 // { (0,phi_i,0) }_{i=0}^{cardPnm1-1} ,
237 // { (0,0,phi_i) }_{i=0}^{cardPnm1-1} ,
238 // { (x,y) . phi_i}_{i=cardPnm2}^{cardPnm1-1}
239 // columns of V1 are expansion of this basis in terms of the orthogonal basis
240 // for P_{n}^3
241
242 // these two loops get the first two sets of basis functions
243 for (ordinal_type i=0;i<cardPnm1;i++) {
244 for (ordinal_type k=0; k<3;k++) {
245 V1(k*cardPn+i,k*cardPnm1+i) = 1.0;
246 }
247 }
248
249 // now I need to integrate { (x,y,z) phi } against the big basis
250 // first, get a cubature rule.
252 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubPoints("Hdiv::Tet::In::cubPoints", myCub.getNumPoints() , spaceDim );
253 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubWeights("Hdiv::Tet::In::cubWeights", myCub.getNumPoints() );
254 myCub.getCubature( cubPoints , cubWeights );
255
256 // tabulate the scalar orthonormal basis at cubature points
257 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtCubPoints("Hdiv::Tet::In::phisAtCubPoints", cardPn , myCub.getNumPoints() );
258 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
259
260 // now do the integration
261 for (ordinal_type i=0;i<dim_PkH;i++) {
262 for (ordinal_type j=0;j<cardPn;j++) { // int (x,y,z) phi_i \cdot (phi_j,0,0)
263 V1(j,cardVecPnm1+i) = 0.0;
264 for (ordinal_type d=0; d< spaceDim; ++d)
265 for (ordinal_type k=0;k<myCub.getNumPoints();k++) {
266 V1(j+d*cardPn,cardVecPnm1+i) +=
267 cubWeights(k) * cubPoints(k,d)
268 * phisAtCubPoints(cardPnm2+i,k)
269 * phisAtCubPoints(j,k);
270 }
271 }
272 }
273
274 // next, apply the RT nodes (rows) to the basis for (P_n)^3 (columns)
275 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
276 V2("Hdiv::Tet::In::V2", card ,cardVecPn);
277
278 const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
279
280 shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
281
282 const int numPtsPerFace = PointTools::getLatticeSize( faceTop ,
283 order+2 ,
284 1 );
285
286 // get the points on the tetrahedron face
287 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> triPts("Hdiv::Tet::In::triPts", numPtsPerFace , 2 );
288
289 // construct lattice
290 const ordinal_type offset = 1;
292 faceTop,
293 order+2,
294 offset,
295 pointType );
296
297 // holds the image of the tet points
298 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> facePts("Hdiv::Tet::In::facePts", numPtsPerFace , spaceDim );
299 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtFacePoints("Hdiv::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace );
300 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceNormal("Hcurl::Tet::In::faceNormal", spaceDim );
301
302 // loop over faces
303 for (ordinal_type face=0;face<numFaces;face++) { // loop over faces
304
305 // these are normal scaled by the appropriate face areas.
307 face ,
308 this->basisCellTopology_ );
309
311 triPts ,
312 2 ,
313 face ,
314 this->basisCellTopology_ );
315
316 // get phi values at face points
317 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints, facePts, order, OPERATOR_VALUE);
318
319 // loop over points (rows of V2)
320 for (ordinal_type j=0;j<numPtsPerFace;j++) {
321
322 const ordinal_type i_card = numPtsPerFace*face+j;
323
324 // loop over orthonormal basis functions (columns of V2)
325 for (ordinal_type k=0;k<cardPn;k++) {
326 // loop over space dimension
327 for (ordinal_type l=0; l<spaceDim; l++)
328 V2(i_card,k+l*cardPn) = faceNormal(l) * phisAtFacePoints(k,j);
329 }
330
331 //save dof coordinates and coefficients
332 for(ordinal_type l=0; l<spaceDim; ++l) {
333 dofCoords(i_card,l) = facePts(j,l);
334 dofCoeffs(i_card,l) = faceNormal(l);
335 }
336
337 tags[i_card][0] = 2; // face dof
338 tags[i_card][1] = face; // face id
339 tags[i_card][2] = j; // local dof id
340 tags[i_card][3] = numPtsPerFace; // total vert dof
341
342 }
343 }
344
345 // remaining nodes point values of each vector component on interior
346 // points of a lattice of degree+2
347 // This way, RT0 --> degree = 1 and internal lattice has no points
348 // RT1 --> degree = 2, and internal lattice has one point (inside of quartic)
349 const ordinal_type numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
350 order + 2 ,
351 1 );
352
353 if (numPtsPerCell > 0) {
354 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
355 internalPoints( "Hdiv::Tet::In::internalPoints", numPtsPerCell , spaceDim );
356 PointTools::getLattice( internalPoints ,
357 this->basisCellTopology_ ,
358 order + 2 ,
359 1 ,
360 pointType );
361
362 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
363 phisAtInternalPoints("Hdiv::Tet::In::phisAtInternalPoints", cardPn , numPtsPerCell );
364 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtInternalPoints , internalPoints , order, OPERATOR_VALUE );
365
366 // copy values into right positions of V2
367 for (ordinal_type j=0;j<numPtsPerCell;j++) {
368
369 const ordinal_type i_card = numFaces*numPtsPerFace+spaceDim*j;
370
371 for (ordinal_type k=0;k<cardPn;k++) {
372 for (ordinal_type l=0;l<spaceDim;l++) {
373 V2(i_card+l,l*cardPn+k) = phisAtInternalPoints(k,j);
374 }
375 }
376
377 //save dof coordinates and coefficients
378 for(ordinal_type d=0; d<spaceDim; ++d) {
379 for(ordinal_type l=0; l<spaceDim; ++l) {
380 dofCoords(i_card+d,l) = internalPoints(j,l);
381 dofCoeffs(i_card+d,l) = (l==d);
382 }
383
384 tags[i_card+d][0] = spaceDim; // elem dof
385 tags[i_card+d][1] = 0; // elem id
386 tags[i_card+d][2] = spaceDim*j+d; // local dof id
387 tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
388 }
389 }
390 }
391
392 // form Vandermonde matrix. Actually, this is the transpose of the VDM,
393 // so we transpose on copy below.
394 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
395 vmat("Hdiv::Tet::In::vmat", card, card),
396 work("Hdiv::Tet::In::work", lwork),
397 ipiv("Hdiv::Tet::In::ipiv", card);
398
399 //vmat' = V2*V1;
400 for(ordinal_type i=0; i< card; ++i) {
401 for(ordinal_type j=0; j< card; ++j) {
402 scalarType s=0;
403 for(ordinal_type k=0; k< cardVecPn; ++k)
404 s += V2(i,k)*V1(k,j);
405 vmat(i,j) = s;
406 }
407 }
408
409 ordinal_type info = 0;
410 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
411
412 lapack.GETRF(card, card,
413 vmat.data(), vmat.stride_1(),
414 (ordinal_type*)ipiv.data(),
415 &info);
416
417 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
418 std::runtime_error ,
419 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRF returns nonzero info." );
420
421 lapack.GETRI(card,
422 vmat.data(), vmat.stride_1(),
423 (ordinal_type*)ipiv.data(),
424 work.data(), lwork,
425 &info);
426
427 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
428 std::runtime_error ,
429 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRI returns nonzero info." );
430
431 for (ordinal_type i=0;i<cardVecPn;++i)
432 for (ordinal_type j=0;j<card;++j){
433 scalarType s=0;
434 for(ordinal_type k=0; k< card; ++k)
435 s += V1(i,k)*vmat(k,j);
436 coeffs(i,j) = s;
437 }
438
439 this->coeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), coeffs);
440 Kokkos::deep_copy(this->coeffs_ , coeffs);
441
442 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
443 Kokkos::deep_copy(this->dofCoords_, dofCoords);
444
445 this->dofCoeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoeffs);
446 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
447
448
449 // set tags
450 {
451 // Basis-dependent initializations
452 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
453 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
454 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
455
456 OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
457
458 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
459 // tags are constructed on host
460 this->setOrdinalTagData(this->tagToOrdinal_,
461 this->ordinalToTag_,
462 tagView,
463 this->basisCardinality_,
464 tagSize,
465 posScDim,
466 posScOrd,
467 posDfOrd);
468 }
469}
470} // namespace Intrepid2
471#endif
Header file for the Intrepid2::CubatureDirectTetDefault class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static void mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties... > refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties... > paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void getReferenceSideNormal(Kokkos::DynRankView< refSideNormalValueType, refSideNormalProperties... > refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
Defines direct integration rules on a tetrahedron.
virtual ordinal_type getNumPoints() const override
Returns the number of cubature points.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
static ordinal_type getLatticeSize(const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties... > points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...