This file is indexed.

/usr/include/trilinos/Intrepid_HGRAD_POLY_C1_FEMDef.hpp is in libtrilinos-intrepid-dev 12.4.2-2.

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
#include "Intrepid_FieldContainer.hpp"
#include "Sacado_Fad_DFad.hpp"


namespace Intrepid{

  template<class Scalar, class ArrayScalar>
  Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::Basis_HGRAD_POLY_C1_FEM(const shards::CellTopology& cellTopology){
    this -> basisCardinality_  = cellTopology.getNodeCount();
    this -> basisDegree_       = 1;
    this -> basisCellTopology_ = cellTopology;
    this -> basisType_         = BASIS_FEM_DEFAULT;
    this -> basisCoordinates_  = COORDINATES_CARTESIAN;
    this -> basisTagsAreSet_   = false;
  }

  
  template<class Scalar, class ArrayScalar>
  void Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::initializeTags(){
    // Basis-dependent initializations
    int tagSize  = 4;        // size of DoF tag, i.e., number of fields in the tag
    int posScDim = 0;        // position in the tag, counting from 0, of the subcell dim 
    int posScOrd = 1;        // position in the tag, counting from 0, of the subcell ordinal
    int posDfOrd = 2;        // position in the tag, counting from 0, of DoF ordinal relative to the subcell
  
    // An array with local DoF tags assigned to the basis functions, in the order of their local enumeration 
    
    int *tags = new int[tagSize * this->getCardinality()];
    for (int i=0;i<this->getCardinality();i++){
      tags[4*i] = 0;
      tags[4*i+1] = i;
      tags[4*i+2] = 0;
      tags[4*i+3] = 1;
    }
    
    
    // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
    Intrepid::setOrdinalTagData(this -> tagToOrdinal_,
				this -> ordinalToTag_,
				tags,
				this -> basisCardinality_,
				tagSize,
				posScDim,
				posScOrd,
				posDfOrd);

    delete [] tags;
  }  
  
  // this is the FEM reference element version, this should not be called for polygonal basis
  // since polygonal basis is defined on physical element.
  template<class Scalar, class ArrayScalar>
  void Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::getValues(ArrayScalar& outputValues,
							       const ArrayScalar& inputPoints,
							       const EOperator operatorType) const{
    TEUCHOS_TEST_FOR_EXCEPTION ( true, std::logic_error,
			 ">>>ERROR (Basis_HGRAD_POLY_C1_FEM): Polygonal basis calling FEM member function");
  }


  template<class Scalar, class ArrayScalar>
  void Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::getValues(ArrayScalar& outputValues,
							       const ArrayScalar& inputPoints,
							       const ArrayScalar& cellVertices,
							       const EOperator operatorType) const{
    
    
    // implement wachspress basis functions
    switch (operatorType) {
    case OPERATOR_VALUE:
      {
	shapeFunctions<Scalar,ArrayScalar>(outputValues,inputPoints,cellVertices);
      }
      break;
    case OPERATOR_GRAD:
      {
	FieldContainer<Sacado::Fad::DFad<Scalar> > dInput(inputPoints.dimension(0),inputPoints.dimension(1));
	for (int i=0;i<inputPoints.dimension(0);i++){
	  for (int j=0;j<2;j++){
	    dInput(i,j) = Sacado::Fad::DFad<Scalar>( inputPoints(i,j));
	    dInput(i,j).diff(j,2);
	  }
	}
	FieldContainer<Sacado::Fad::DFad<Scalar> > cellVerticesFad(cellVertices.dimension(0),cellVertices.dimension(1));
	for (int i=0;i<cellVertices.dimension(0);i++){
	  for (int j=0;j<cellVertices.dimension(1);j++){
	    cellVerticesFad(i,j) = Sacado::Fad::DFad<Scalar>( cellVertices(i,j) );
	  }
	}
	
	FieldContainer<Sacado::Fad::DFad<Scalar> > dOutput(outputValues.dimension(0),outputValues.dimension(1));
	
	shapeFunctions<Sacado::Fad::DFad<Scalar>, FieldContainer<Sacado::Fad::DFad<Scalar> > >(dOutput,dInput,cellVerticesFad);
	
	for (int i=0;i<outputValues.dimension(0);i++){
	  for (int j=0;j<outputValues.dimension(1);j++){
	    for (int k=0;k<outputValues.dimension(2);k++){
	      outputValues(i,j,k) = dOutput(i,j).dx(k);
	    }
	  }
	}
      }
      break;
    case OPERATOR_D1:
    case OPERATOR_D2:
    case OPERATOR_D3:
    case OPERATOR_D4:
    case OPERATOR_D5:
    case OPERATOR_D6:
    case OPERATOR_D7:
    case OPERATOR_D8:
    case OPERATOR_D9:
    case OPERATOR_D10:
      {
	TEUCHOS_TEST_FOR_EXCEPTION ( true, std::invalid_argument, 
			     ">>> ERROR (Basis_HGRAD_POLY_C1_FEM): operator not implemented yet");
      }
      break;
    default:
      TEUCHOS_TEST_FOR_EXCEPTION( !( Intrepid::isValidOperator(operatorType)), std::invalid_argument,
			  ">>> ERROR (Basis_HGRAD_POLY_C1_FEM): Invalid operator type");
      break;
    }

  }

  
  template<class Scalar, class ArrayScalar>
  template<class Scalar1, class ArrayScalar1>
  void Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::shapeFunctions(ArrayScalar1& outputValues,
								    const ArrayScalar1& inputPoints,
								    const ArrayScalar1& cellVertices)const{
    int numPoints = inputPoints.dimension(0);
    FieldContainer<Scalar1> weightFunctions( this->basisCardinality_,numPoints );
    evaluateWeightFunctions<Scalar1, ArrayScalar1>(weightFunctions,inputPoints,cellVertices);
    for (int pointIndex = 0;pointIndex<outputValues.dimension(1);pointIndex++){
      Scalar1 denominator=0;
      
      for (int k=0;k<weightFunctions.dimension(0);k++){
	denominator += weightFunctions(k,pointIndex);
      }
      for (int k=0;k<outputValues.dimension(0);k++){
	outputValues(k,pointIndex) = weightFunctions(k,pointIndex)/denominator;
      }
    }
  }
							    


  template<class Scalar, class ArrayScalar>
  template<class Scalar1, class ArrayScalar1>
  Scalar1 Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::computeArea(const ArrayScalar1& p1,
								    const ArrayScalar1& p2,
								    const ArrayScalar1& p3) const{
    Scalar1 area = 0.5*(p1(0)*(p2(1)-p3(1))
			-p2(0)*(p1(1)-p3(1))
			+p3(0)*(p1(1)-p2(1)));
    return area;
  }
  

  template<class Scalar, class ArrayScalar>
  template<class Scalar1, class ArrayScalar1>
  void Basis_HGRAD_POLY_C1_FEM<Scalar, ArrayScalar>::evaluateWeightFunctions(ArrayScalar1& outputValues,
									     const ArrayScalar1& inputValues,
									     const ArrayScalar1& cellVertices)const{
    
    
    int spaceDim = this->basisCellTopology_.getDimension();
    for (int k=0;k<outputValues.dimension(0);k++){
      
      // compute a_k for each weight function
      // need a_k = A(p_i,p_j,p_k) where nodes i and j are adjacent to node k
      int adjIndex1 = -1, adjIndex2 = -1;
      for (int i = 0;i < this->basisCellTopology_.getEdgeCount();i++){
	if ( this->basisCellTopology_.getNodeMap(1,i,0) == k )
	  adjIndex1 = this->basisCellTopology_.getNodeMap(1,i,1);
	else if ( this->basisCellTopology_.getNodeMap(1,i,1) == k )
	  adjIndex2 = this->basisCellTopology_.getNodeMap(1,i,0);
      }
      TEUCHOS_TEST_FOR_EXCEPTION( (adjIndex1 == -1 || adjIndex2 == -1), std::invalid_argument,
			  ">>> ERROR (Intrepid_HGRAD_POLY_C1_FEM): cannot find adjacent nodes when evaluating Wachspress weight function.");
      FieldContainer<Scalar1> p1(spaceDim);
      FieldContainer<Scalar1> p2(spaceDim);
      FieldContainer<Scalar1> p3(spaceDim);
      for (int i=0;i<spaceDim;i++){
	p1(i) = cellVertices(adjIndex1,i);
	p2(i) = cellVertices(k,i);
	p3(i) = cellVertices(adjIndex2,i);
      }
      Scalar1 a_k = computeArea<Scalar1, ArrayScalar1>(p1,p2,p3);
      // now compute prod_{ij!=k} r_ij
      for (int pointIndex = 0;pointIndex < inputValues.dimension(0);pointIndex++){
	Scalar1 product = a_k;
	for (int edgeIndex = 0;edgeIndex < this->basisCellTopology_.getEdgeCount();edgeIndex++){
	  int edgeNode1 = this->basisCellTopology_.getNodeMap(1,edgeIndex,0);
	  int edgeNode2 = this->basisCellTopology_.getNodeMap(1,edgeIndex,1);
	  if ( edgeNode1 != k && edgeNode2 != k ){
	    for (int i=0;i<spaceDim;i++){
	      p1(i) = inputValues(pointIndex,i);
	      p2(i) = cellVertices(edgeNode1,i);
	      p3(i) = cellVertices(edgeNode2,i);
	    }
	    product *= computeArea<Scalar1, ArrayScalar1>(p1,p2,p3);
	  }
	}
	outputValues(k,pointIndex) = product;
      }
    }
  }
} // namespace Intrepid