File: CubatureTensorTests.cpp

package info (click to toggle)
trilinos 13.2.0-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 837,948 kB
  • sloc: cpp: 3,466,824; ansic: 433,701; fortran: 168,838; python: 102,712; xml: 66,353; sh: 38,380; makefile: 25,246; perl: 7,516; f90: 6,358; csh: 4,161; objc: 2,620; lex: 1,484; lisp: 810; javascript: 552; yacc: 515; awk: 364; ml: 281; php: 145
file content (191 lines) | stat: -rw-r--r-- 7,413 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
// @HEADER
// ************************************************************************
//
//                           Intrepid2 Package
//                 Copyright (2007) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Kyungjoo Kim  (kyukim@sandia.gov),
//                    Mauro Perego  (mperego@sandia.gov), or
//                    Nate Roberts  (nvrober@sandia.gov)
//
// ************************************************************************
// @HEADER

/** \file   BasisValuesTests.cpp
    \brief  Tests to verify that the version of Basis::getValues() that takes the new BasisValues and TensorPoints arguments produces the same results as the one that takes raw Kokkos DynRankViews as arguments.
    \author Created by N.V. Roberts.
 */

#include "Teuchos_UnitTestHarness.hpp"

#include "Intrepid2_DefaultCubatureFactory.hpp"
#include "Intrepid2_ScalarView.hpp"
#include "Intrepid2_Types.hpp"
#include "Intrepid2_TestUtils.hpp"

#include <Kokkos_Core.hpp>

using namespace Intrepid2;

namespace
{
  using namespace Intrepid2;

  void testTensorPointCubature(const shards::CellTopology &cellTopo, const int &quadratureDegree,
                               const double &relTol, const double &absTol, Teuchos::FancyOStream &out, bool &success)
  {
    using DeviceType = DefaultTestDeviceType;
    using PointScalar = double;
    using WeightScalar = double;
    using CubatureType   = Cubature<DeviceType,PointScalar,WeightScalar>;
    using PointViewType  = typename CubatureType::PointViewTypeAllocatable;
    using WeightViewType = typename CubatureType::WeightViewTypeAllocatable;
    
    DefaultCubatureFactory cub_factory;
    auto cellTopoKey = cellTopo.getKey();
    auto quadrature = cub_factory.create<DeviceType, PointScalar, WeightScalar>(cellTopoKey, quadratureDegree);
    ordinal_type numRefPoints = quadrature->getNumPoints();
    const int spaceDim = cellTopo.getDimension();
    
    PointViewType points("quadrature points ref cell", numRefPoints, spaceDim);
    WeightViewType weights("quadrature weights ref cell", numRefPoints);
    
    quadrature->getCubature(points, weights);
    
    TensorPoints<PointScalar,DeviceType> tensorPoints;
    TensorData<WeightScalar,DeviceType>  tensorWeights;
    
    using CubatureTensorType = CubatureTensor<DeviceType,PointScalar,WeightScalar>;
    CubatureTensorType* tensorQuadrature = dynamic_cast<CubatureTensorType*>(quadrature.get());

    if (tensorQuadrature)
    {
      tensorPoints  = tensorQuadrature->allocateCubaturePoints();
      tensorWeights = tensorQuadrature->allocateCubatureWeights();
      tensorQuadrature->getCubature(tensorPoints, tensorWeights);
    }
    else
    {
      std::vector<ViewType<PointScalar,DeviceType>> pointComponents {points};
      tensorPoints = TensorPoints<PointScalar,DeviceType>(pointComponents);
      Data<WeightScalar,DeviceType> weightData(weights);
      std::vector<Data<WeightScalar,DeviceType>> weightComponents {weightData};
      tensorWeights = TensorData<WeightScalar,DeviceType>(weightComponents);
    }
    
    printView(points, out, "Points being tested");
    printFunctor2(tensorPoints, out, "tensorPoints");
    
    // points and tensorPoints should be identical, no roundoff: use 0.0 tolerances
    testFloatingEquality2(points,  tensorPoints,    0.0,    0.0, out, success, "points",  "tensorPoints");
    testFloatingEquality1(weights,tensorWeights, relTol, absTol, out, success, "weights", "tensorWeights");
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Line )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Line<> >());
    
    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Quadrilateral )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<> >());
    
    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Hexahedron )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<> >());

    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Tetrahedron )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<> >());
    
    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Triangle )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<> >());
    
    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }

  TEUCHOS_UNIT_TEST( CubatureTensor, Wedge )
  {
    shards::CellTopology cellTopo = shards::CellTopology(shards::getCellTopologyData<shards::Wedge<> >());
    
    const double relTol=1e-13;
    const double absTol=1e-13;
    
    for (int polyOrder=1; polyOrder<5; polyOrder++)
    {
      testTensorPointCubature(cellTopo, polyOrder, relTol, absTol, out, success);
    }
  }
} // namespace