This file is indexed.

/usr/share/pyshared/FIAT/raviart_thomas.py is in python-fiat 1.0.0-1.

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
import expansions, polynomial_set, quadrature, reference_element, dual_set, \
    quadrature, finite_element, functional
import numpy

def RTSpace( ref_el , deg ):
    """Constructs a basis for the the Raviart-Thomas space
    (P_k)^d + P_k x"""
    sd = ref_el.get_spatial_dimension()

    vec_Pkp1 = polynomial_set.ONPolynomialSet( ref_el , deg+1 , (sd,) )

    dimPkp1 = expansions.polynomial_dimension( ref_el , deg+1 )
    dimPk = expansions.polynomial_dimension( ref_el , deg )
    dimPkm1 = expansions.polynomial_dimension( ref_el , deg-1 )

    vec_Pk_indices = reduce( lambda a,b: a+b , \
                             [ range(i*dimPkp1,i*dimPkp1+dimPk) \
                               for i in range(sd) ] )
    vec_Pk_from_Pkp1 = vec_Pkp1.take( vec_Pk_indices )

    Pkp1 = polynomial_set.ONPolynomialSet( ref_el , deg + 1 )
    PkH = Pkp1.take( range(dimPkm1,dimPk) )

    Q = quadrature.make_quadrature( ref_el , 2 * deg + 2 )

    # have to work on this through "tabulate" interface
    # first, tabulate PkH at quadrature points
    Qpts = numpy.array( Q.get_points() )
    Qwts = numpy.array( Q.get_weights() )

    zero_index = tuple( [ 0 for i in range(sd) ] )

    PkH_at_Qpts = PkH.tabulate( Qpts )[zero_index]
    Pkp1_at_Qpts = Pkp1.tabulate( Qpts )[zero_index]

    PkHx_coeffs = numpy.zeros( (PkH.get_num_members() , \
                                sd, \
                                Pkp1.get_num_members()) , "d" )

    import time
    t1 = time.time()
    for i in range( PkH.get_num_members() ):
        for j in range( sd ):
            fooij = PkH_at_Qpts[i,:] * Qpts[:,j] * Qwts
            PkHx_coeffs[i,j,:] = numpy.dot( Pkp1_at_Qpts , fooij )

    PkHx = polynomial_set.PolynomialSet( ref_el , \
                                         deg , \
                                         deg + 1 , \
                                         vec_Pkp1.get_expansion_set() , \
                                         PkHx_coeffs , \
                                         vec_Pkp1.get_dmats() )

    return polynomial_set.polynomial_set_union_normalized( vec_Pk_from_Pkp1 , PkHx )

class RTDualSet( dual_set.DualSet ):
    """Dual basis for Raviart-Thomas elements consisting of point
    evaluation of normals on facets of codimension 1 and internal
    moments against polynomials"""
    def __init__( self , ref_el , degree ):
        entity_ids = {}
        nodes = []

        sd = ref_el.get_spatial_dimension()


        t = ref_el.get_topology()

        # codimension 1 facets

        for i in range( len( t[sd-1] ) ):
            pts_cur = ref_el.make_points( sd - 1 , i , sd + degree )
            for j in range( len( pts_cur ) ):
                pt_cur = pts_cur[j]
                f = functional.PointScaledNormalEvaluation( ref_el , i , \
                                                            pt_cur )
                nodes.append( f )

        # internal nodes.  Let's just use points at a lattice
        if degree > 0:
            cpe = functional.ComponentPointEvaluation
            pts = ref_el.make_points( sd , 0 , degree + sd )
            for d in range( sd ):
                for i in range( len( pts ) ):
                    l_cur = cpe( ref_el , d , (sd,) , pts[i] )
                    nodes.append( l_cur )


#            Q = quadrature.make_quadrature( ref_el , 2 * ( degree + 1 ) )
#            qpts = Q.get_points()
#            Pkm1 = polynomial_set.ONPolynomialSet( ref_el , degree - 1 )
#            zero_index = tuple( [ 0 for i in range( sd ) ]  )
#            Pkm1_at_qpts = Pkm1.tabulate( qpts )[ zero_index ]

#            for d in range( sd ):
#                for i in range( Pkm1_at_qpts.shape[0] ):
#                    phi_cur = Pkm1_at_qpts[i,:]
#                    l_cur = functional.IntegralMoment( ref_el , Q , \
#                                                       phi_cur , (d,) , (sd,) )
#                    nodes.append( l_cur )


        entity_ids = {}

        # sets vertices (and in 3d, edges) to have no nodes
        for i in range( sd - 1 ):
            entity_ids[i] = {}
            for j in range( len( t[i] ) ):
                entity_ids[i][j] = []

        cur = 0

        # set codimension 1 (edges 2d, faces 3d) dof
        pts_facet_0 = ref_el.make_points( sd - 1 , 0 , sd + degree )
        pts_per_facet = len( pts_facet_0 )

        entity_ids[sd-1] = {}
        for i in range( len( t[sd-1] ) ):
            entity_ids[sd-1][i] = range( cur , cur + pts_per_facet )
            cur += pts_per_facet

        # internal nodes, if applicable
        entity_ids[sd] = {0: []}

        if degree > 0:
            num_internal_nodes = expansions.polynomial_dimension( ref_el , \
                                                                  degree - 1 )
            entity_ids[sd][0] = range( cur , cur + num_internal_nodes * sd )

        dual_set.DualSet.__init__( self , nodes , ref_el , entity_ids )


class RaviartThomas( finite_element.FiniteElement ):
    """The Raviart-Thomas finite element"""
    def __init__( self , ref_el , q ):

        degree = q - 1
        poly_set = RTSpace( ref_el , degree )
        dual = RTDualSet( ref_el , degree )
        finite_element.FiniteElement.__init__( self , poly_set , dual , degree,
                                               mapping="contravariant piola")
        return


if __name__=="__main__":
    T = reference_element.UFCTriangle()
    sd = T.get_spatial_dimension()

    for k in range(6):
        RT = RaviartThomas( T , k )


#    RTfs = RT.get_nodal_basis()

#    pts = T.make_lattice( 1 )
#    print pts

#    zero_index = tuple( [ 0 for i in range(sd) ] )
#
#    RTvals = RTfs.tabulate( pts )[zero_index]

#    print RTvals