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 "Panzer_PureBasis.hpp"
51 
52 #include "Teuchos_FancyOStream.hpp"
53 
54 template<typename EvalT,typename Traits>
57  const Teuchos::ParameterList& p)
58 {
59  dof_name = (p.get< std::string >("DOF Name"));
60 
61  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
62  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
63  else
64  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
65 
66  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
67  Teuchos::RCP<PHX::DataLayout> vector_layout = basis->functional_grad;
68 
69  // some sanity checks
70  TEUCHOS_ASSERT(basis->isVectorBasis());
71 
72  result = PHX::MDField<ScalarT,Cell,BASIS>(dof_name,basis_layout);
73  this->addEvaluatedField(result);
74 
75  vector_values = PHX::MDField<ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector",vector_layout);
76  tangents = PHX::MDField<ScalarT,Cell,BASIS,Dim>(dof_name+"_Tangents",vector_layout);
77  this->addDependentField(vector_values);
78  this->addDependentField(tangents);
79 
80  this->setName("Project To Edges");
81 }
82 
83 // **********************************************************************
84 template<typename EvalT,typename Traits>
88 {
89  // setup the field data object
90  this->utils.setFieldData(result,fm);
91  this->utils.setFieldData(vector_values,fm);
92  this->utils.setFieldData(tangents,fm);
93 
94  num_pts = vector_values.dimension(1);
95  num_dim = vector_values.dimension(2);
96 
97  TEUCHOS_ASSERT(vector_values.dimension(1) == tangents.dimension(1));
98  TEUCHOS_ASSERT(vector_values.dimension(2) == tangents.dimension(2));
99 }
100 
101 // **********************************************************************
102 template<typename EvalT,typename Traits>
105 {
106  const shards::CellTopology & parentCell = *basis->getCellTopology();
107  const int intDegree = basis->order();
108  TEUCHOS_ASSERT(intDegree == 1);
109  Intrepid2::DefaultCubatureFactory<double,Kokkos::DynRankView<double,PHX::Device>,Kokkos::DynRankView<double,PHX::Device>> quadFactory;
110  Teuchos::RCP< Intrepid2::Cubature<double,Kokkos::DynRankView<double,PHX::Device>,Kokkos::DynRankView<double,PHX::Device>> > edgeQuad;
111 
112  // Collect the reference edge information. For now, do nothing with the quadPts.
113  const unsigned num_edges = parentCell.getEdgeCount();
114  std::vector<double> refEdgeWt(num_edges, 0.0);
115  for (unsigned e=0; e<num_edges; e++) {
116  edgeQuad = quadFactory.create(parentCell.getCellTopologyData(1,e), intDegree);
117  const int numQPoints = edgeQuad->getNumPoints();
118  Kokkos::DynRankView<double,PHX::Device> quadWts("quadWts",numQPoints);
119  Kokkos::DynRankView<double,PHX::Device> quadPts("quadPts",numQPoints,num_dim);
120  edgeQuad->getCubature(quadPts,quadWts);
121  for (int q=0; q<numQPoints; q++)
122  refEdgeWt[e] += quadWts(q);
123  }
124 
125  // Loop over the edges of the workset cells.
126  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
127  for (int p = 0; p < num_pts; ++p) {
128  result(cell,p) = ScalarT(0.0);
129  for (int dim = 0; dim < num_dim; ++dim)
130  result(cell,p) += vector_values(cell,p,dim) * tangents(cell,p,dim);
131  result(cell,p) *= refEdgeWt[p];
132  }
133  }
134 
135 }
136 
137 #endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.