File: DirectSumBasisTests.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 (121 lines) | stat: -rw-r--r-- 5,154 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
// @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   DirectSumBasisTests.cpp
    \brief  Tests against the DirectSumBasis class.
    \author Created by N.V. Roberts.
 */

#include "Teuchos_UnitTestHarness.hpp"

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

#include <Kokkos_Core.hpp>

namespace
{
  using namespace Intrepid2;
// #pragma mark DirectSumBasis group
  TEUCHOS_UNIT_TEST( DirectSumBasis, AllocateBasisValues )
  {
    using Basis = HierarchicalBasisFamily<DefaultTestDeviceType>::HDIV_HEX;
    
    Basis basis(1,1,1);
    
    int quadratureDegree = 1;
    using DeviceType = DefaultTestDeviceType;
    using PointScalar = double;
    using WeightScalar = double;
    DefaultCubatureFactory cub_factory;
    using CubatureType   = Cubature<DeviceType,PointScalar,WeightScalar>;
    using PointViewType  = typename CubatureType::PointViewTypeAllocatable;
    using WeightViewType = typename CubatureType::WeightViewTypeAllocatable;
    auto cellTopo = basis.getBaseCellTopology();
    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);
    }
    
    // test OPERATOR_DIV
    auto basisValues = basis.allocateBasisValues(tensorPoints, OPERATOR_DIV);
    
    const int basisCardinality = basis.getCardinality();
    TEST_EQUALITY(basisCardinality, basisValues.extent_int(0));
    TEST_EQUALITY(numRefPoints, basisValues.extent_int(1));
    TEST_EQUALITY(2, basisValues.rank());
    
    int family1Fields = basisValues.numFieldsInFamily(0);
    int family2Fields = basisValues.numFieldsInFamily(1);
    int family3Fields = basisValues.numFieldsInFamily(2);
    
    TEST_EQUALITY(basisCardinality, family1Fields + family2Fields + family3Fields);
    
    // test OPERATOR_VALUE
    basisValues = basis.allocateBasisValues(tensorPoints, OPERATOR_VALUE);
    TEST_EQUALITY(basisCardinality, basisValues.extent_int(0));
    TEST_EQUALITY(numRefPoints, basisValues.extent_int(1));
    TEST_EQUALITY(spaceDim, basisValues.extent_int(2));
    TEST_EQUALITY(3, basisValues.rank());
  }
} // namespace