Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_GatherFields_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_GATHER_FIELDS_IMPL_HPP
12#define PANZER_STK_GATHER_FIELDS_IMPL_HPP
13
14#include "Teuchos_Assert.hpp"
15#include "Phalanx_DataLayout.hpp"
16
18
19#include "Teuchos_FancyOStream.hpp"
20
21// **********************************************************************
22// Specialization: Residual
23// **********************************************************************
24
25template<typename EvalT, typename Traits>
27 GatherFields(const Teuchos::RCP<const STK_Interface> & mesh,const Teuchos::ParameterList& p)
28{
29 using panzer::Cell;
30 using panzer::NODE;
31
32 mesh_ = mesh;
33
34 const std::vector<std::string>& names =
35 *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
36
37 Teuchos::RCP<panzer::BasisIRLayout> basis =
38 p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis");
39 isConstant_ = basis->getBasis()->getElementSpace()==panzer::PureBasis::CONST;
40
41 gatherFields_.resize(names.size());
42 stkFields_.resize(names.size());
43 for (std::size_t fd = 0; fd < names.size(); ++fd) {
44 gatherFields_[fd] =
45 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
46 this->addEvaluatedField(gatherFields_[fd]);
47 }
48
49 this->setName("Gather STK Fields");
50}
51
52// **********************************************************************
53template<typename EvalT, typename Traits>
55postRegistrationSetup(typename Traits::SetupData /* d */,
57{
58 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
59 std::string fieldName = gatherFields_[fd].fieldTag().name();
60
61 stkFields_[fd] = mesh_->getMetaData()->get_field<double>(stk::topology::NODE_RANK, fieldName);
62
63 if(stkFields_[fd]==0) {
64 std::stringstream ss;
65 mesh_->printMetaData(ss);
66 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
67 "panzer_stk::GatherFields: STK field " << "\"" << fieldName << "\" "
68 "not found.\n STK meta data follows: \n\n" << ss.str());
69 }
70 }
71}
72
73// **********************************************************************
74template<typename EvalT, typename Traits>
76evaluateFields(typename Traits::EvalData workset)
77{
78 const std::vector<stk::mesh::Entity> & localElements = *mesh_->getElementsOrderedByLID();
79
80 // for convenience pull out some objects from workset
81 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
82
83 // gather operation for each cell in workset
84 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
85 std::size_t cellLocalId = localCellIds[worksetCellIndex];
86 stk::mesh::Entity const* relations = mesh_->getBulkData()->begin_nodes(localElements[cellLocalId]);
87
88 // loop over the fields to be gathered
89 for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
90 VariableField * field = stkFields_[fieldIndex];
91
92 std::size_t basisCnt = gatherFields_[fieldIndex].extent(1);
93
94 if(isConstant_) {
95 // loop over basis functions and fill the fields
96 (gatherFields_[fieldIndex])(worksetCellIndex,0) = *stk::mesh::field_data(*field, localElements[cellLocalId]);
97 }
98 else {
99 // loop over basis functions and fill the fields
100 for(std::size_t basis=0;basis<basisCnt;basis++) {
101 stk::mesh::Entity node = relations[basis];
102 (gatherFields_[fieldIndex])(worksetCellIndex,basis) = *stk::mesh::field_data(*field, node);
103 }
104 }
105 }
106 }
107}
108
109#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.
void evaluateFields(typename Traits::EvalData d)
panzer_stk::STK_Interface::SolutionFieldType VariableField
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)