Panzer
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
adapters-stk
src
evaluators
Panzer_STK_ScatterCellAvgVector_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_CELL_AVG_VECTOR_IMPL_HPP
44
#define PANZER_STK_SCATTER_CELL_AVG_VECTOR_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_IntegrationRule.hpp
"
55
#include "
Panzer_CommonArrayFactories.hpp
"
56
57
#include "Teuchos_FancyOStream.hpp"
58
#include "Teuchos_ArrayRCP.hpp"
59
60
namespace
panzer_stk
{
61
62
template
<
typename
EvalT,
typename
Traits>
63
ScatterCellAvgVector<EvalT, Traits>::
64
ScatterCellAvgVector
(
65
const
Teuchos::ParameterList& p) :
66
mesh_
(p.get<
Teuchos
::RCP<
STK_Interface
> >(
"Mesh"
))
67
{
68
using
panzer::Cell
;
69
using
panzer::Point
;
70
using
panzer::Dim
;
71
72
std::string scatterName = p.get<std::string>(
"Scatter Name"
);
73
74
if
(p.isParameter(
"Variable Scale Factors Map"
))
75
varScaleFactors_
= p.get<Teuchos::RCP<std::map<std::string,double>>>(
"Variable Scale Factors Map"
);
76
77
const
std::vector<std::string> & names =
78
*(p.get< Teuchos::RCP< std::vector<std::string> > >(
"Field Names"
));
79
80
Teuchos::RCP<panzer::IntegrationRule> intRule =
81
p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR"
);
82
83
// build dependent fields
84
scatterFields_
.resize(names.size());
85
stkFields_
.resize(names.size());
86
for (std::size_t fd = 0; fd < names.size(); ++fd)
87
{
88
scatterFields_
[fd] = PHX::MDField<const ScalarT,Cell,Point,Dim>(names[fd],intRule->dl_vector);
89
this->addDependentField(
scatterFields_
[fd]);
90
}
91
92
// setup a dummy field to evaluate
93
PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(
new
PHX::MDALayout<panzer::Dummy>(0)));
94
this->addEvaluatedField(scatterHolder);
95
96
this->setName(scatterName+
": STK-Scatter Cell Vectors"
);
97
}
98
99
100
template
<
typename
EvalT,
typename
Traits>
101
void
102
ScatterCellAvgVector<EvalT, Traits>::
103
postRegistrationSetup
(
104
typename
Traits::SetupData
/* d */
,
105
PHX::FieldManager<Traits>
&
/* fm */
)
106
{
107
for
(std::size_t fd = 0; fd <
scatterFields_
.size(); ++fd)
108
{
109
std::string fieldName =
scatterFields_
[fd].fieldTag().name();
110
111
stkFields_
[fd] =
mesh_
->getMetaData()->get_field<
VariableField
>(stk::topology::ELEMENT_RANK, fieldName);
112
}
113
}
114
115
116
template
<
typename
EvalT,
typename
Traits>
117
void
118
ScatterCellAvgVector<EvalT, Traits>::
119
evaluateFields
(
120
typename
Traits::EvalData workset)
121
{
122
panzer::MDFieldArrayFactory
af(
""
,
true
);
123
124
// for convenience pull out some objects from workset
125
const
std::vector<std::size_t> & localCellIds = this->
wda
(workset).cell_local_ids;
126
std::string blockId = this->
wda
(workset).block_id;
127
std::string d_mod[3] = {
"X"
,
"Y"
,
"Z"
};
128
129
// loop over the number of vector fields requested for exodus output
130
for
(std::size_t fieldIndex = 0; fieldIndex <
scatterFields_
.size(); fieldIndex++)
131
{
132
PHX::MDField<const ScalarT,panzer::Cell,panzer::Point,panzer::Dim> & field =
scatterFields_
[fieldIndex];
133
std::string fieldName = field.fieldTag().name();
134
int
numCells = field.extent(0);
135
int
numPoints = field.extent(1);
136
int
numDims = field.extent(2);
137
138
for
(
int
dim = 0; dim < numDims; dim++)
139
{
140
// std::vector<double> average(numCells,0.0);
141
PHX::MDField<double,panzer::Cell,panzer::NODE> average = af.
buildStaticArray
<double,
panzer::Cell
,
panzer::NODE
>(
""
,numCells,1);
142
143
// write to double field
144
for
(
int
i = 0; i < numCells; i++)
// loop over cells
145
{
146
average(i,0) = 0.0;
147
for
(
int
j = 0; j < numPoints; j++)
// loop over IPs
148
average(i,0) += Sacado::scalarValue(field(i,j,dim));
149
150
average(i,0) /= numPoints;
151
}
152
153
double
scalef = 1.0;
154
155
if
(!
varScaleFactors_
.is_null())
156
{
157
std::map<std::string,double> *tmp_sfs =
varScaleFactors_
.get();
158
if
(tmp_sfs->find(fieldName) != tmp_sfs->end())
159
scalef = (*tmp_sfs)[fieldName];
160
}
161
162
mesh_
->setCellFieldData(fieldName+d_mod[dim],blockId,localCellIds,average.get_view(),scalef);
163
164
}
165
}
166
167
}
168
169
}
// end panzer_stk
170
171
#endif
Panzer_CommonArrayFactories.hpp
Panzer_IntegrationRule.hpp
PHX::FieldManager
Definition
Panzer_BCStrategy_Base.hpp:53
panzer::EvaluatorWithBaseImpl< Traits >::wda
WorksetDetailsAccessor wda
Definition
Panzer_Evaluator_WithBaseImpl.hpp:63
panzer::MDFieldArrayFactory
Definition
Panzer_CommonArrayFactories.hpp:86
panzer::MDFieldArrayFactory::buildStaticArray
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
Definition
Panzer_CommonArrayFactories_impl.hpp:167
panzer_stk::STK_Interface
Definition
Panzer_STK_Interface.hpp:104
panzer_stk::ScatterCellAvgVector::varScaleFactors_
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_
Definition
Panzer_STK_ScatterCellAvgVector_decl.hpp:106
panzer_stk::ScatterCellAvgVector::scatterFields_
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point, panzer::Dim > > scatterFields_
Definition
Panzer_STK_ScatterCellAvgVector_decl.hpp:108
panzer_stk::ScatterCellAvgVector::mesh_
Teuchos::RCP< STK_Interface > mesh_
Definition
Panzer_STK_ScatterCellAvgVector_decl.hpp:109
panzer_stk::ScatterCellAvgVector::postRegistrationSetup
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Definition
Panzer_STK_ScatterCellAvgVector_impl.hpp:103
panzer_stk::ScatterCellAvgVector::ScatterCellAvgVector
ScatterCellAvgVector(const Teuchos::ParameterList &p)
Definition
Panzer_STK_ScatterCellAvgVector_impl.hpp:64
panzer_stk::ScatterCellAvgVector::evaluateFields
void evaluateFields(typename Traits::EvalData d)
Definition
Panzer_STK_ScatterCellAvgVector_impl.hpp:119
panzer_stk::ScatterCellAvgVector::VariableField
panzer_stk::STK_Interface::VectorFieldType VariableField
Definition
Panzer_STK_ScatterCellAvgVector_decl.hpp:98
panzer_stk::ScatterCellAvgVector::stkFields_
std::vector< VariableField * > stkFields_
Definition
Panzer_STK_ScatterCellAvgVector_decl.hpp:111
Teuchos
panzer_stk
Definition
Panzer_STK_GatherFields_decl.hpp:58
panzer::NODE
BASIS NODE
Definition
Panzer_Dimension.hpp:77
panzer::Cell
Definition
Panzer_Dimension.hpp:79
panzer::Dim
Definition
Panzer_Dimension.hpp:74
panzer::Point
Definition
Panzer_Dimension.hpp:78
Generated by
1.17.0