Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_TensorToStdVector_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_TENSOR_TO_STD_VECTOR_IMPL_HPP
44#define PANZER_TENSOR_TO_STD_VECTOR_IMPL_HPP
45
46#include <string>
47
48namespace panzer {
49
50//**********************************************************************
51template<typename EvalT, typename Traits>
54 const Teuchos::ParameterList& p)
55{
56 Teuchos::RCP<PHX::DataLayout> vector_dl =
57 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
58
59 Teuchos::RCP<PHX::DataLayout> tensor_dl =
60 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Tensor");
61
62 const std::vector<std::string>& vector_names =
63 *(p.get< Teuchos::RCP<const std::vector<std::string> > >("Vector Names"));
64
65 vector_fields.resize(vector_names.size());
66 for (std::size_t i=0; i < vector_names.size(); ++i)
67 vector_fields[i] =
68 PHX::MDField<ScalarT,Cell,Point,Dim>(vector_names[i], vector_dl);
69
70 tensor_field =
71 PHX::MDField<const ScalarT,Cell,Point,Dim,Dim>(p.get<std::string>
72 ("Tensor Name"), tensor_dl);
73
74 this->addDependentField(tensor_field);
75
76 for (std::size_t i=0; i < vector_fields.size(); ++i)
77 this->addEvaluatedField(vector_fields[i]);
78
79 std::string n = "TensorToStdVector: " + tensor_field.fieldTag().name();
80 this->setName(n);
81}
82
83//**********************************************************************
84template<typename EvalT, typename Traits>
85void
88 typename Traits::EvalData workset)
89{
90
91 typedef typename PHX::MDField<ScalarT,Cell,Point,Dim>::size_type size_type;
92
93 // Loop over cells
94 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
95
96 // Loop over points
97 for (size_type pt = 0; pt < tensor_field.extent(1); ++pt) {
98
99 // Loop over vectors
100 for (std::size_t vec = 0; vec < vector_fields.size(); ++vec) {
101
102 // Loop over spatial dimensions
103 for (std::size_t dim = 0; dim < tensor_field.extent(2); ++dim) {
104
105 vector_fields[vec](cell,pt,dim) = tensor_field(cell,pt,vec,dim);
106
107 }
108 }
109 }
110 }
111}
112
113//**********************************************************************
114
115}
116
117#endif
void evaluateFields(typename Traits::EvalData d)
TensorToStdVector(const Teuchos::ParameterList &p)
int num_cells
DEPRECATED - use: numCells()