Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_ScatterCellAvgQuantity_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_CELL_AVG_QUANTITY_IMPL_HPP
12#define PANZER_STK_SCATTER_CELL_AVG_QUANTITY_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
24
25#include "Teuchos_FancyOStream.hpp"
26#include "Teuchos_ArrayRCP.hpp"
27
28namespace panzer_stk {
29
30template<typename EvalT, typename Traits>
33 const Teuchos::ParameterList& p) :
34 mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
35{
36 using panzer::Cell;
37 using panzer::Point;
38
39 std::string scatterName = p.get<std::string>("Scatter Name");
40
41 // if it's there pull the ouptut scaling
42 if (p.isParameter("Variable Scale Factors Map"))
43 varScaleFactors_ = p.get<Teuchos::RCP<std::map<std::string,double>>>("Variable Scale Factors Map");
44
45 const std::vector<std::string> & names =
46 *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
47
48 Teuchos::RCP<panzer::IntegrationRule> intRule =
49 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
50
51 // build dependent fields
52 scatterFields_.resize(names.size());
53 stkFields_.resize(names.size());
54 for (std::size_t fd = 0; fd < names.size(); ++fd) {
55 scatterFields_[fd] =
56 PHX::MDField<const ScalarT,Cell,Point>(names[fd],intRule->dl_scalar);
57 this->addDependentField(scatterFields_[fd]);
58 }
59
60 // setup a dummy field to evaluate
61 PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
62 this->addEvaluatedField(scatterHolder);
63
64 this->setName(scatterName+": STK-Scatter Cell Fields");
65}
66
67template<typename EvalT, typename Traits>
68void
71 typename Traits::SetupData /* d */,
73{
74 for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd) {
75 std::string fieldName = scatterFields_[fd].fieldTag().name();
76
77 stkFields_[fd] = mesh_->getMetaData()->get_field<double>(stk::topology::ELEMENT_RANK, fieldName);
78 }
79}
80
81template<typename EvalT, typename Traits>
82void
85 typename Traits::EvalData workset)
86{
88
89 // for convenience pull out some objects from workset
90 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
91 std::string blockId = this->wda(workset).block_id;
92
93 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
94 auto field = scatterFields_[fieldIndex].get_static_view();
95 auto average = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.extent(0),1).get_static_view();
96 // write to double field
97 Kokkos::parallel_for("ScatterCellAvgQuantity",field.extent(0), KOKKOS_LAMBDA(int i) {
98 for(unsigned j=0; j<field.extent(1);j++)
99 average(i,0) += Sacado::scalarValue(field(i,j));
100 average(i,0) /= field.extent(1);
101 });
102 Kokkos::fence();
103
104 std::string varname = scatterFields_[fieldIndex].fieldTag().name();
105 double scalef = 1.0;
106
107 if (!varScaleFactors_.is_null())
108 {
109 std::map<std::string,double> *tmp_sfs = varScaleFactors_.get();
110 if(tmp_sfs->find(varname) != tmp_sfs->end())
111 scalef = (*tmp_sfs)[varname];
112 }
113
114 mesh_->setCellFieldData(varname,blockId,localCellIds,average,scalef);
115 }
116}
117
118} // end panzer_stk
119
120#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
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point > > scatterFields_