Panzer  Version of the Day
Panzer_ProjectToEdges_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_PROJECT_TO_EDGES_IMPL_HPP
44 #define PANZER_PROJECT_TO_EDGES_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Intrepid2_Cubature.hpp"
50 #include "Intrepid2_DefaultCubatureFactory.hpp"
51 #include "Intrepid2_FunctionSpaceTools.hpp"
52 #include "Intrepid2_OrientationTools.hpp"
53 
54 #include "Panzer_PureBasis.hpp"
56 #include "Kokkos_ViewFactory.hpp"
57 
58 #include "Teuchos_FancyOStream.hpp"
59 
60 template<typename EvalT,typename Traits>
63  const Teuchos::ParameterList& p)
64 {
65  dof_name = (p.get< std::string >("DOF Name"));
66 
67  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
68  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
69  else
70  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
71 
72  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
73  Teuchos::RCP<PHX::DataLayout> vector_layout = basis->functional_grad;
74 
75  // some sanity checks
76  TEUCHOS_ASSERT(basis->isVectorBasis());
77 
78  result = PHX::MDField<ScalarT,Cell,BASIS>(dof_name,basis_layout);
79  this->addEvaluatedField(result);
80 
81  tangents = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Tangents",vector_layout);
82  this->addDependentField(tangents);
83 
84  vector_values.resize(1);
85  vector_values[0] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector",vector_layout);
86  this->addDependentField(vector_values[0]);
87 
88  this->setName("Project To Edges");
89 }
90 
91 // **********************************************************************
92 template<typename EvalT,typename Traits>
95  PHX::FieldManager<Traits>& /* fm */)
96 {
97  num_edges = vector_values[0].extent(1);
98  num_dim = vector_values[0].extent(2);
99 
100  TEUCHOS_ASSERT(vector_values[0].extent(1) == tangents.extent(1));
101  TEUCHOS_ASSERT(vector_values[0].extent(2) == tangents.extent(2));
102 }
103 
104 // **********************************************************************
105 template<typename EvalT,typename Traits>
108 {
109  // Restricting HCurl field, multiplied by the tangent to the edge, into HVol on the edges.
110  // This code assumes affine mapping and the projection into 1 quadrature point for each edge,
111  // which is identified with the edge. This makes sense only for low order bases, for which
112  // HVol is constant
113 
114  //TODO: make this work w/ high order basis
115  const int intDegree = basis->order();
116  TEUCHOS_ASSERT(intDegree == 1);
117 
118  //TODO: use parallel for over cells
119  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
120  for (int p = 0; p < num_edges; ++p) {
121  result(cell,p) = ScalarT(0.0);
122  for (int dim = 0; dim < num_dim; ++dim)
123  result(cell,p) += vector_values[0](cell,p,dim) * tangents(cell,p,dim);
124  }
125  }
126 }
127 
128 #endif
int num_cells
DEPRECATED - use: numCells()
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we&#39;re performing.
void evaluateFields(typename Traits::EvalData d)