Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_SubcellSum_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_SubcellSum_impl_hpp__
12#define __Panzer_SubcellSum_impl_hpp__
13
14#include "Panzer_PureBasis.hpp"
17
18#include "Phalanx_DataLayout_MDALayout.hpp"
19
20namespace panzer {
21
22//**********************************************************************
23template<typename EvalT, typename Traits>
26 const Teuchos::ParameterList& p)
27 : evaluateOnClosure_(false)
28{
29 Teuchos::RCP<Teuchos::ParameterList> valid_params = this->getValidParameters();
30 p.validateParameters(*valid_params);
31
32 const std::string inName = p.get<std::string>("Field Name");
33 const std::string outName = p.get<std::string>("Sum Name");
34 Teuchos::RCP<const PureBasis> basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
35 multiplier = p.get<double>("Multiplier");
36 if(p.isType<bool>("Evaluate On Closure"))
37 evaluateOnClosure_ = p.get<bool>("Evaluate On Closure");
38
39 inField = PHX::MDField<const ScalarT,Cell,BASIS>( inName, basis->functional);
40 outField = PHX::MDField<ScalarT,Cell>( outName, basis->cell_data);
41
42 this->addDependentField(inField);
43 this->addEvaluatedField(outField);
44
45 // build a field pattern object so that looking up closure indices is easy
46 fieldPattern_ = Teuchos::rcp(new Intrepid2FieldPattern(basis->getIntrepid2Basis<PHX::exec_space,double,double>()));
47
48 std::string n = "SubcellSum: " + outField.fieldTag().name();
49 this->setName(n);
50}
51
52//**********************************************************************
53template<typename EvalT, typename Traits>
54void
57 typename Traits::EvalData workset)
58{
59 std::vector<int> indices;
60
61 // figure out which indices to sum (this can be made more efficient by
62 // simply saving the indices and only updating if the subcell dimension
63 // and index changes)
64 if(evaluateOnClosure_)
65 fieldPattern_->getSubcellClosureIndices(workset.subcell_dim,this->wda(workset).subcell_index,indices);
66 else
67 indices = fieldPattern_->getSubcellIndices(workset.subcell_dim,this->wda(workset).subcell_index);
68
69 auto outField_h = Kokkos::create_mirror_view(outField.get_static_view());
70 auto inField_h = Kokkos::create_mirror_view(inField.get_static_view());
71 Kokkos::deep_copy(inField_h, inField.get_static_view());
72 for(index_t c=0;c<workset.num_cells;c++) {
73 outField_h(c) = 0.0; // initialize field
74
75 // sum over all relevant indices for this subcell
76 for(std::size_t i=0;i<indices.size();i++)
77 outField_h(c) += inField_h(c,indices[i]);
78
79 // scale by what ever the user wants
80 outField_h(c) *= multiplier;
81 }
82 Kokkos::deep_copy(outField.get_static_view(), outField_h);
83}
84
85//**********************************************************************
86template<typename EvalT, typename TRAITS>
87Teuchos::RCP<Teuchos::ParameterList>
89{
90 Teuchos::RCP<Teuchos::ParameterList> p = Teuchos::rcp(new Teuchos::ParameterList);
91 p->set<std::string>("Sum Name", "?");
92 p->set<std::string>("Field Name", "?");
93 p->set<double>("Multiplier",1.0);
94 p->set<bool>("Evaluate On Closure",false);
95
96 Teuchos::RCP<const panzer::PureBasis> basis;
97 p->set("Basis", basis);
98
99 return p;
100}
101
102//**********************************************************************
103
104}
105
106#endif
double multiplier
The scalar multiplier out in front of the integral ( ).
SubcellSum(const Teuchos::ParameterList &p)
PHX::MDField< ScalarT, Cell > outField
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const panzer::FieldPattern > fieldPattern_
void evaluateFields(typename Traits::EvalData d)
PHX::MDField< const ScalarT, Cell, BASIS > inField
int subcell_dim
DEPRECATED - use: getSubcellDimension()
int num_cells
DEPRECATED - use: numCells()