This file is indexed.

/usr/include/trilinos/Intrepid_CubatureControlVolumeBoundaryDef.hpp is in libtrilinos-intrepid-dev 12.12.1-5.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
// @HEADER
// ************************************************************************
//
//                           Intrepid 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 Pavel Bochev  (pbboche@sandia.gov)
//                    Denis Ridzal  (dridzal@sandia.gov), or
//                    Kara Peterson (kjpeter@sandia.gov)
//
// ************************************************************************
// @HEADER

/** \file   Intrepid_CubatureControlVolumeBoundaryDef.hpp
    \brief  Definition file for the Intrepid::CubatureControlVolumeBoundary class.
    \author Created by K. Peterson, P. Bochev and D. Ridzal.
*/

#ifndef INTREPID_CUBATURE_CONTROLVOLUMEBOUNDARYDEF_HPP
#define INTREPID_CUBATURE_CONTROLVOLUMEBOUNDARYDEF_HPP

namespace Intrepid{

template<class Scalar, class ArrayPoint, class ArrayWeight>
CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::CubatureControlVolumeBoundary(const Teuchos::RCP<const shards::CellTopology> & cellTopology, int cellSide)
{
  // topology of primary cell 
  primaryCellTopo_ = cellTopology;

  // get topology of sub-control volume (will be quad or hex depending on dimension)
  const CellTopologyData &myCellData =
  (primaryCellTopo_->getDimension() > 2) ? *shards::getCellTopologyData<shards::Hexahedron<8> >() :
                                           *shards::getCellTopologyData<shards::Quadrilateral<4> >();

  subCVCellTopo_ = Teuchos::rcp(new shards::CellTopology(&myCellData));

  degree_ = 1;

  cubDimension_ = primaryCellTopo_->getDimension();

  sideIndex_ = cellSide;

  // one control volume boundary cubature point per subcell node (for now)
  numPoints_ = primaryCellTopo_->getNodeCount(primaryCellTopo_->getDimension()-1,cellSide);

}

template<class Scalar, class ArrayPoint, class ArrayWeight>
void CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::getCubature(ArrayPoint& cubPoints,
		                                                               ArrayWeight& cubWeights) const
{
    TEUCHOS_TEST_FOR_EXCEPTION( (true), std::logic_error,
                      ">>> ERROR (CubatureControlVolumeBoundary): Cubature defined in physical space calling method for reference space cubature.");
}

template<class Scalar, class ArrayPoint, class ArrayWeight>
void CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::getCubature(ArrayPoint& cubPoints,
  		                                                               ArrayWeight& cubWeights,
                                                                               ArrayPoint& cellCoords) const
{
  // get array dimensions
  int numCells         = cellCoords.dimension(0);
  int numNodesPerCell  = cellCoords.dimension(1);
  int spaceDim         = cellCoords.dimension(2);
  int numNodesPerSubCV = subCVCellTopo_->getNodeCount();

  // get sub-control volume coordinates (one sub-control volume per node of primary cell)
  Intrepid::FieldContainer<Scalar> subCVCoords(numCells,numNodesPerCell,numNodesPerSubCV,spaceDim);
  Intrepid::CellTools<Scalar>::getSubCVCoords(subCVCoords,cellCoords,*(primaryCellTopo_));

  // define subcontrol volume side index corresponding to primary cell side
  int numPrimarySideNodes = primaryCellTopo_->getNodeCount(spaceDim-1,sideIndex_);  
  int numCVSideNodes = subCVCellTopo_->getNodeCount(spaceDim-1,0);  
  int numPrimarySides = primaryCellTopo_->getSubcellCount(spaceDim-1);
  Intrepid::FieldContainer<int> CVSideonBoundary(numPrimarySides,numCVSideNodes);

  switch(primaryCellTopo_->getKey() ) {

      case shards::Triangle<3>::key:
      case shards::Quadrilateral<4>::key:

         for (int iside=0; iside<numPrimarySides; iside++) {
              CVSideonBoundary(iside,0) = 0; CVSideonBoundary(iside,1) = 3;
          }

      break;

      case shards::Hexahedron<8>::key:

         // sides 0-3
         for (int iside=0; iside<4; iside++) {
              CVSideonBoundary(iside,0) = 0; CVSideonBoundary(iside,1) = 3;
              CVSideonBoundary(iside,2) = 3; CVSideonBoundary(iside,3) = 0;
         }
         // side 4
         CVSideonBoundary(4,0) = 4; CVSideonBoundary(4,1) = 4;
         CVSideonBoundary(4,2) = 4; CVSideonBoundary(4,3) = 4;

         // side 5
         CVSideonBoundary(5,0) = 5; CVSideonBoundary(5,1) = 5;
         CVSideonBoundary(5,2) = 5; CVSideonBoundary(5,3) = 5;

      break;

      case shards::Tetrahedron<4>::key:

         CVSideonBoundary(0,0) = 0; CVSideonBoundary(0,1) = 3; CVSideonBoundary(0,2) = 0;
         CVSideonBoundary(1,0) = 0; CVSideonBoundary(1,1) = 3; CVSideonBoundary(1,2) = 3;
         CVSideonBoundary(2,0) = 3; CVSideonBoundary(2,1) = 4; CVSideonBoundary(2,2) = 0;
         CVSideonBoundary(3,0) = 4; CVSideonBoundary(3,1) = 4; CVSideonBoundary(3,2) = 4;

      break;

      default:
        TEUCHOS_TEST_FOR_EXCEPTION( true, std::invalid_argument,
                            ">>> ERROR (CubatureControlVolumeBoundary: invalid cell topology.");
    } // cell key

    Intrepid::FieldContainer<double> sideCenterLocal(1,spaceDim-1);
    for (int idim = 0; idim < spaceDim-1; idim++){
        sideCenterLocal(0,idim) = 0.0;
    }

    // get side cubature points
    for (int icell = 0; icell < numCells; icell++)
    {

      for (int inode=0; inode<numPrimarySideNodes; inode++){

         int cvind = primaryCellTopo_->getNodeMap(spaceDim-1,sideIndex_,inode);
         int cvside = CVSideonBoundary(sideIndex_,inode);

         Intrepid::FieldContainer<double> cubpoint(spaceDim);
         for (int idim=0; idim<spaceDim; idim++) { 
            for (int icvnode=0; icvnode<numCVSideNodes; icvnode++) { 
               int cvnode = subCVCellTopo_->getNodeMap(spaceDim-1,cvside,icvnode);
               cubpoint(idim) = cubpoint(idim) + subCVCoords(icell,cvind,cvnode,idim);       
            }
            cubPoints(icell,inode,idim) = cubpoint(idim)/numCVSideNodes; 
         }

         // map side center to reference subcell
          Intrepid::FieldContainer<Scalar> refSidePoints(1,spaceDim);
          Intrepid::CellTools<Scalar>::mapToReferenceSubcell(refSidePoints,
                                        sideCenterLocal,
                                        spaceDim-1, cvside, *(subCVCellTopo_));

         // array of sub-control volume coordinates
          Intrepid::FieldContainer<Scalar> cellCVCoords(1, numNodesPerSubCV, spaceDim);
          for (int icvnode = 0; icvnode < numNodesPerSubCV; icvnode++){
              for (int idim = 0; idim < spaceDim; idim++){
                   cellCVCoords(0,icvnode,idim) = subCVCoords(icell,cvind,icvnode,idim);
              }
          }

         // calculate Jacobian at side centers
          Intrepid::FieldContainer<Scalar> subCVsideJacobian(1, 1, spaceDim, spaceDim);
          Intrepid::FieldContainer<Scalar> subCVsideJacobianDet(1, 1);
          Intrepid::CellTools<Scalar>::setJacobian(subCVsideJacobian, refSidePoints, cellCVCoords, *(subCVCellTopo_));
          Intrepid::CellTools<Scalar>::setJacobianDet(subCVsideJacobianDet, subCVsideJacobian);

         // calculate Jacobian at side centers
          Intrepid::FieldContainer<Scalar> measure(1, 1);
          Intrepid::FieldContainer<Scalar> weights(1, 1);
          if (spaceDim == 3){
             weights(0,0) = 4.0;
             Intrepid::FunctionSpaceTools::computeFaceMeasure<Scalar>(measure,subCVsideJacobian,weights,cvside,*(subCVCellTopo_));
          }
          else if (spaceDim == 2){
             weights(0,0) = 2.0;
             Intrepid::FunctionSpaceTools::computeEdgeMeasure<Scalar>(measure,subCVsideJacobian,weights,cvside,*(subCVCellTopo_));
          }

          cubWeights(icell,inode) = measure(0,0);

       } // end loop over primary side nodes

 } // end cell loop

} // end getCubature
    

template<class Scalar, class ArrayPoint, class ArrayWeight>
int CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::getNumPoints()const{
  return numPoints_;
} // end getNumPoints

template<class Scalar, class ArrayPoint, class ArrayWeight>
int CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::getDimension()const{
  return cubDimension_;
} // end getDimension

template<class Scalar, class ArrayPoint, class ArrayWeight>
void CubatureControlVolumeBoundary<Scalar,ArrayPoint,ArrayWeight>::getAccuracy(std::vector<int>& accuracy)const{
  accuracy.assign(1,degree_);
} // end getAccuracy

} // end namespace Intrepid

#endif