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// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_STK_SCATTER_VECTOR_FIELDS_IMPL_HPP
12#define PANZER_STK_SCATTER_VECTOR_FIELDS_IMPL_HPP
13
14#include "Teuchos_Assert.hpp"
15
16#include "Phalanx_config.hpp"
17#include "Phalanx_Evaluator_Macros.hpp"
18#include "Phalanx_MDField.hpp"
19#include "Phalanx_DataLayout.hpp"
20#include "Phalanx_DataLayout_MDALayout.hpp"
21
22#include "Panzer_Traits.hpp"
24
25#include "Teuchos_FancyOStream.hpp"
26
27namespace panzer_stk {
28
29template<typename EvalT, typename Traits>
32 const Teuchos::ParameterList& p) :
33 mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
34{
35 TEUCHOS_ASSERT(false);
36}
37
38template <typename EvalT,typename TraitsT>
40ScatterVectorFields(const std::string & scatterName,
41 const Teuchos::RCP<STK_Interface> mesh,
42 const Teuchos::RCP<const panzer::PointRule> & pointRule,
43 const std::vector<std::string> & names,
44 const std::vector<double> & scaling)
45 : mesh_(mesh)
46 , scaling_(scaling)
47{
48 using panzer::Cell;
49 using panzer::IP;
50 using panzer::Dim;
51
52 spatialDimension_ = pointRule->spatial_dimension;
53
54 // this evaluator assumes you are evaluating at the cell centroid only
55 TEUCHOS_ASSERT(pointRule->num_points==1);
56
57 // build dependent fields
58 names_ = names;
59 scatterFields_.resize(names.size());
60 for (std::size_t fd = 0; fd < names.size(); ++fd) {
61 scatterFields_[fd] =
62 PHX::MDField<const ScalarT,Cell,IP,Dim>(names_[fd]+"_"+pointRule->getName(),pointRule->dl_vector);
63 this->addDependentField(scatterFields_[fd]);
64 }
65
66 // setup a dummy field to evaluate
67 PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
68 this->addEvaluatedField(scatterHolder);
69
70 this->setName(scatterName+": STK-Scatter Vector Fields");
71}
72
73template<typename EvalT, typename Traits>
74void
77 typename Traits::EvalData /* workset */)
78{
79 TEUCHOS_ASSERT(false);
80}
81
82template < >
85{
87
88 std::vector<std::string> dimStrings(3);
89 dimStrings[0] = "X";
90 dimStrings[1] = "Y";
91 dimStrings[2] = "Z";
92
93 // for convenience pull out some objects from workset
94 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
95 std::string blockId = this->wda(workset).block_id;
96
97 for(int d=0;d<spatialDimension_;d++) {
98 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
99 PHX::MDField<const ScalarT,panzer::Cell,panzer::IP,panzer::Dim> & field = scatterFields_[fieldIndex];
100 std::string fieldName = names_[fieldIndex]+dimStrings[d];
101
102 PHX::MDField<double,panzer::Cell,panzer::NODE> cellValue
103 = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.extent(0),1);
104
105 // scaline field value only if the scaling parameter is specified, otherwise use 1.0
106 double scaling = (scaling_.size()>0) ? scaling_[fieldIndex] : 1.0;
107
108 auto cellValue_v = cellValue.get_static_view();
109 auto field_v = field.get_static_view();
110 Kokkos::parallel_for(field_v.extent(0), KOKKOS_LAMBDA (int i) {
111 cellValue_v(i,0) = field_v(i,0,d);
112 });
113 Kokkos::fence();
114
115 // add in vector value at d^th dimension
116 mesh_->setCellFieldData(fieldName,blockId,localCellIds,cellValue.get_view(),scaling);
117 }
118 }
119}
120
121} // end panzer_stk
122
123#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_