File: test_util.hpp

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 (132 lines) | stat: -rw-r--r-- 5,102 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
// @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), or
//                    Mauro Perego  (mperego@sandia.gov)
//
// ************************************************************************
// @HEADER


/** \file
    \brief  Test utils (compute volume, monomial, integral)
    \author Kyungjoo Kim
*/

#ifndef __INTREPID2_TEST_UTILS_HPP__
#define __INTREPID2_TEST_UTILS_HPP__

#include "Intrepid2_config.h"

#include "Intrepid2_Types.hpp"
#include "Intrepid2_Utils.hpp"

namespace Intrepid2 {

  namespace Test {

    template<typename cubWeightViewType>
    typename cubWeightViewType::value_type
    computeRefVolume(const ordinal_type      numPoints,
                     const cubWeightViewType cubWeights) {
      typename cubWeightViewType::value_type r_val = 0.0;
      Kokkos::fence();
      auto cubWeights_host = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), cubWeights);
      for (auto i=0;i<numPoints;++i) {
        r_val += cubWeights_host(i);
      }

      return r_val;
    }

    // Monomial evaluation.
    // in 1D, for point p(x)    : x^xDeg
    // in 2D, for point p(x,y)  : x^xDeg * y^yDeg
    // in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
    template<typename ValueType, 
             typename PointViewType>
    ValueType computeMonomial(PointViewType p, 
                              const ordinal_type xDeg, 
                              const ordinal_type yDeg = 0, 
                              const ordinal_type zDeg = 0) {
      ValueType r_val = 1.0;
      const ordinal_type polydeg[3] = { xDeg, yDeg, zDeg };
      const auto dim = p.extent(0);
      Kokkos::fence();
      
      Kokkos::DynRankView<typename PointViewType::non_const_value_type,
                          typename PointViewType::execution_space> p_device("p_device", p.extent(0));
      Kokkos::deep_copy(p_device, p);
      auto p_host = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), p_device);
      for (size_type i=0;i<dim;++i) 
        r_val *= std::pow(p_host(i),polydeg[i]);
      
      return r_val;
    }

    // Computes integrals of monomials 
    template<typename ValueType,
             typename cubatureType,
             typename cubPointViewType,
             typename cubWeightViewType>
    ValueType computeIntegralOfMonomial(cubatureType cub,
                                        cubPointViewType cubPoints,
                                        cubWeightViewType cubWeights,
                                        const ordinal_type xDeg,
                                        const ordinal_type yDeg = 0,
                                        const ordinal_type zDeg = 0) {
      ValueType r_val = 0.0;

      // get cubature 
      cub.getCubature(cubPoints, cubWeights);
      auto cubWeights_host = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), cubWeights);
      Kokkos::fence();

      const auto dim  = cub.getDimension();
      const auto npts = cub.getNumPoints();
      typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;

      for (auto i=0;i<npts;++i) {
        const auto pt = Kokkos::subdynrankview(cubPoints, i, range_type(0, dim));
        r_val += computeMonomial<ValueType>(pt, xDeg, yDeg, zDeg)*cubWeights_host(i);
      }

      return r_val;
    }
  }
}

#endif