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//
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_GATHER_FIELDS_IMPL_HPP
44#define PANZER_STK_GATHER_FIELDS_IMPL_HPP
45
46#include "Teuchos_Assert.hpp"
47#include "Phalanx_DataLayout.hpp"
48
50
51#include "Teuchos_FancyOStream.hpp"
52
53// **********************************************************************
54// Specialization: Residual
55// **********************************************************************
56
57template<typename EvalT, typename Traits>
59 GatherFields(const Teuchos::RCP<const STK_Interface> & mesh,const Teuchos::ParameterList& p)
60{
61 using panzer::Cell;
62 using panzer::NODE;
63
64 mesh_ = mesh;
65
66 const std::vector<std::string>& names =
67 *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
68
69 Teuchos::RCP<panzer::BasisIRLayout> basis =
70 p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis");
71 isConstant_ = basis->getBasis()->getElementSpace()==panzer::PureBasis::CONST;
72
73 gatherFields_.resize(names.size());
74 stkFields_.resize(names.size());
75 for (std::size_t fd = 0; fd < names.size(); ++fd) {
76 gatherFields_[fd] =
77 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
78 this->addEvaluatedField(gatherFields_[fd]);
79 }
80
81 this->setName("Gather STK Fields");
82}
83
84// **********************************************************************
85template<typename EvalT, typename Traits>
87postRegistrationSetup(typename Traits::SetupData /* d */,
89{
90 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
91 std::string fieldName = gatherFields_[fd].fieldTag().name();
92
93 stkFields_[fd] = mesh_->getMetaData()->get_field<VariableField>(stk::topology::NODE_RANK, fieldName);
94
95 if(stkFields_[fd]==0) {
96 std::stringstream ss;
97 mesh_->printMetaData(ss);
98 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
99 "panzer_stk::GatherFields: STK field " << "\"" << fieldName << "\" "
100 "not found.\n STK meta data follows: \n\n" << ss.str());
101 }
102 }
103}
104
105// **********************************************************************
106template<typename EvalT, typename Traits>
108evaluateFields(typename Traits::EvalData workset)
109{
110 const std::vector<stk::mesh::Entity> & localElements = *mesh_->getElementsOrderedByLID();
111
112 // for convenience pull out some objects from workset
113 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
114
115 // gather operation for each cell in workset
116 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
117 std::size_t cellLocalId = localCellIds[worksetCellIndex];
118 stk::mesh::Entity const* relations = mesh_->getBulkData()->begin_nodes(localElements[cellLocalId]);
119
120 // loop over the fields to be gathered
121 for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
122 VariableField * field = stkFields_[fieldIndex];
123
124 std::size_t basisCnt = gatherFields_[fieldIndex].extent(1);
125
126 if(isConstant_) {
127 // loop over basis functions and fill the fields
128 (gatherFields_[fieldIndex])(worksetCellIndex,0) = *stk::mesh::field_data(*field, localElements[cellLocalId]);
129 }
130 else {
131 // loop over basis functions and fill the fields
132 for(std::size_t basis=0;basis<basisCnt;basis++) {
133 stk::mesh::Entity node = relations[basis];
134 (gatherFields_[fieldIndex])(worksetCellIndex,basis) = *stk::mesh::field_data(*field, node);
135 }
136 }
137 }
138 }
139}
140
141#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)