Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_ScatterFields_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_FIELDS_IMPL_HPP
12#define PANZER_STK_SCATTER_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
23#include "Panzer_Traits.hpp"
24
25#include "Teuchos_FancyOStream.hpp"
26
27namespace panzer_stk {
28
29template <typename EvalT,typename TraitsT>
31ScatterFields(const std::string & scatterName,
32 const Teuchos::RCP<STK_Interface> mesh,
33 const Teuchos::RCP<const panzer::PureBasis> & basis,
34 const std::vector<std::string> & names)
35{
36 std::vector<double> scaling; // empty
37
38 initialize(scatterName,mesh,basis,names,scaling);
39}
40
41template <typename EvalT,typename TraitsT>
43ScatterFields(const std::string & scatterName,
44 const Teuchos::RCP<STK_Interface> mesh,
45 const Teuchos::RCP<const panzer::PureBasis> & basis,
46 const std::vector<std::string> & names,
47 const std::vector<double> & scaling)
48{
49 initialize(scatterName,mesh,basis,names,scaling);
50}
51
52template <typename EvalT,typename TraitsT>
54initialize(const std::string & scatterName,
55 const Teuchos::RCP<STK_Interface> mesh,
56 const Teuchos::RCP<const panzer::PureBasis> & basis,
57 const std::vector<std::string> & names,
58 const std::vector<double> & scaling)
59{
60 using panzer::Cell;
61 using panzer::NODE;
62
63 mesh_ = mesh;
64 scaling_ = scaling;
65
66 bool correctScaling = (names.size()==scaling.size()) || (scaling.size()==0);
67 TEUCHOS_TEST_FOR_EXCEPTION(!correctScaling,std::invalid_argument,
68 "panzer_stk::ScatterFields evaluator requites a consistent number of scaling parameters (equal to the number of field names) "
69 "or an empty \"Field Scaling\" vector");
70
71 // build dependent fields
72 scatterFields_.resize(names.size());
73 for (std::size_t fd = 0; fd < names.size(); ++fd) {
74 scatterFields_[fd] =
75 PHX::MDField<const ScalarT,Cell,NODE>(names[fd],basis->functional);
76 this->addDependentField(scatterFields_[fd]);
77 }
78
79 // determine if this is a cell field or not
80 cellFields_ = basis->getElementSpace()==panzer::PureBasis::CONST;
81
82 // setup a dummy field to evaluate
83 PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
84 this->addEvaluatedField(scatterHolder);
85
86 this->setName(scatterName+": STK-Scatter Fields");
87}
88
89template <typename EvalT,typename TraitsT>
91postRegistrationSetup(typename TraitsT::SetupData /* d */,
93{
94 for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
95 std::string fieldName = scatterFields_[fd].fieldTag().name();
96}
97
98template <typename EvalT,typename TraitsT>
100evaluateFields(typename TraitsT::EvalData /* d */)
101{
102 TEUCHOS_ASSERT(false);
103}
104
105template < >
108{
109 // for convenience pull out some objects from workset
110 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
111 std::string blockId = this->wda(workset).block_id;
112
113 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
114 // scaline field value only if the scaling parameter is specified, otherwise use 1.0
115 double scaling = (scaling_.size()>0) ? scaling_[fieldIndex] : 1.0;
116
117 // write field to the STK mesh object
118 if(!cellFields_)
119 mesh_->setSolutionFieldData(scatterFields_[fieldIndex].fieldTag().name(),blockId,
120 localCellIds,scatterFields_[fieldIndex].get_view(),scaling);
121 else
122 mesh_->setCellFieldData(scatterFields_[fieldIndex].fieldTag().name(),blockId,
123 localCellIds,scatterFields_[fieldIndex].get_view(),scaling);
124 }
125}
126
127} // end panzer_stk
128
129#endif
void initialize(const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PureBasis > &basis, const std::vector< std::string > &names, const std::vector< double > &scaling)
void evaluateFields(typename TraitsT::EvalData d)
ScatterFields(const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PureBasis > &basis, const std::vector< std::string > &names)
void postRegistrationSetup(typename TraitsT::SetupData d, PHX::FieldManager< TraitsT > &fm)