Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_ScatterVectorFields_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_STK_SCATTER_VECTOR_FIELDS_IMPL_HPP
44#define PANZER_STK_SCATTER_VECTOR_FIELDS_IMPL_HPP
45
46#include "Teuchos_Assert.hpp"
47
48#include "Phalanx_config.hpp"
49#include "Phalanx_Evaluator_Macros.hpp"
50#include "Phalanx_MDField.hpp"
51#include "Phalanx_DataLayout.hpp"
52#include "Phalanx_DataLayout_MDALayout.hpp"
53
54#include "Panzer_Traits.hpp"
56
57#include "Teuchos_FancyOStream.hpp"
58
59namespace panzer_stk {
60
61template<typename EvalT, typename Traits>
64 const Teuchos::ParameterList& p) :
65 mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
66{
67 TEUCHOS_ASSERT(false);
68}
69
70template <typename EvalT,typename TraitsT>
72ScatterVectorFields(const std::string & scatterName,
73 const Teuchos::RCP<STK_Interface> mesh,
74 const Teuchos::RCP<const panzer::PointRule> & pointRule,
75 const std::vector<std::string> & names,
76 const std::vector<double> & scaling)
77 : mesh_(mesh)
78 , scaling_(scaling)
79{
80 using panzer::Cell;
81 using panzer::IP;
82 using panzer::Dim;
83
84 spatialDimension_ = pointRule->spatial_dimension;
85
86 // this evaluator assumes you are evaluating at the cell centroid only
87 TEUCHOS_ASSERT(pointRule->num_points==1);
88
89 // build dependent fields
90 names_ = names;
91 scatterFields_.resize(names.size());
92 for (std::size_t fd = 0; fd < names.size(); ++fd) {
93 scatterFields_[fd] =
94 PHX::MDField<const ScalarT,Cell,IP,Dim>(names_[fd]+"_"+pointRule->getName(),pointRule->dl_vector);
95 this->addDependentField(scatterFields_[fd]);
96 }
97
98 // setup a dummy field to evaluate
99 PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
100 this->addEvaluatedField(scatterHolder);
101
102 this->setName(scatterName+": STK-Scatter Vector Fields");
103}
104
105template<typename EvalT, typename Traits>
106void
109 typename Traits::EvalData /* workset */)
110{
111 TEUCHOS_ASSERT(false);
112}
113
114template < >
117{
118 panzer::MDFieldArrayFactory af("",true);
119
120 std::vector<std::string> dimStrings(3);
121 dimStrings[0] = "X";
122 dimStrings[1] = "Y";
123 dimStrings[2] = "Z";
124
125 // for convenience pull out some objects from workset
126 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
127 std::string blockId = this->wda(workset).block_id;
128
129 for(int d=0;d<spatialDimension_;d++) {
130 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
131 PHX::MDField<const ScalarT,panzer::Cell,panzer::IP,panzer::Dim> & field = scatterFields_[fieldIndex];
132 std::string fieldName = names_[fieldIndex]+dimStrings[d];
133
134 PHX::MDField<double,panzer::Cell,panzer::NODE> cellValue
135 = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.extent(0),1);
136
137 // scaline field value only if the scaling parameter is specified, otherwise use 1.0
138 double scaling = (scaling_.size()>0) ? scaling_[fieldIndex] : 1.0;
139
140 auto cellValue_v = cellValue.get_static_view();
141 auto field_v = field.get_static_view();
142 Kokkos::parallel_for(field_v.extent(0), KOKKOS_LAMBDA (int i) {
143 cellValue_v(i,0) = field_v(i,0,d);
144 });
145 Kokkos::fence();
146
147 // add in vector value at d^th dimension
148 mesh_->setCellFieldData(fieldName,blockId,localCellIds,cellValue.get_view(),scaling);
149 }
150 }
151}
152
153} // end panzer_stk
154
155#endif
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
ScatterVectorFields(const Teuchos::ParameterList &p)
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > > scatterFields_
void evaluateFields(typename Traits::EvalData d)