Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_CrossProduct_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_EVALUATOR_CrossProduct_IMPL_HPP
12#define PANZER_EVALUATOR_CrossProduct_IMPL_HPP
13
14#include <string>
15
16#include "Panzer_PointRule.hpp"
18
19#include "Teuchos_RCP.hpp"
20
21namespace panzer {
22
23//**********************************************************************
24template<typename EvalT, typename Traits>
27 const Teuchos::ParameterList& p)
28{
29 std::string result_name = p.get<std::string>("Result Name");
30 std::string vec_a_name = p.get<std::string>("Vector A Name");
31 std::string vec_b_name = p.get<std::string>("Vector B Name");
32
33 const Teuchos::RCP<const panzer::PointRule> pr =
34 p.get< Teuchos::RCP<const panzer::PointRule> >("Point Rule");
35
36 // use a scalar field only if dimension is 2D
37 useScalarField = (pr->spatial_dimension==2);
38
39 if(!useScalarField)
40 vec_a_cross_vec_b = PHX::MDField<ScalarT>(result_name, pr->dl_vector);
41 else
42 vec_a_cross_vec_b = PHX::MDField<ScalarT>(result_name, pr->dl_scalar);
43
44 vec_a = PHX::MDField<const ScalarT>(vec_a_name, pr->dl_vector);
45 vec_b = PHX::MDField<const ScalarT>(vec_b_name, pr->dl_vector);
46
47 this->addEvaluatedField(vec_a_cross_vec_b);
48 this->addDependentField(vec_a);
49 this->addDependentField(vec_b);
50
51 std::string n = "CrossProduct: " + result_name + " = " + vec_a_name + " . " + vec_b_name;
52 this->setName(n);
53}
54
55//**********************************************************************
56template<typename EvalT, typename Traits>
57void
60 typename Traits::SetupData /* sd */,
62{
63 num_pts = vec_a.extent(1);
64 num_dim = vec_a.extent(2);
65
66 TEUCHOS_ASSERT(vec_a.extent(1) == vec_b.extent(1));
67 TEUCHOS_ASSERT(vec_a.extent(2) == vec_b.extent(2));
68}
69
70//**********************************************************************
71template<typename EvalT, typename Traits>
72void
75 typename Traits::EvalData workset)
76{
77 if(useScalarField) {
78 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
79 for (int p = 0; p < num_pts; ++p) {
80 vec_a_cross_vec_b(cell,p) = vec_a(cell,p,0)*vec_b(cell,p,1)-vec_a(cell,p,1)*vec_b(cell,p,0);
81 }
82 }
83 }
84 else {
85 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
86 for (int p = 0; p < num_pts; ++p) {
87 vec_a_cross_vec_b(cell,p,0) = vec_a(cell,p,1)*vec_b(cell,p,2)-vec_a(cell,p,2)*vec_b(cell,p,1);
88 vec_a_cross_vec_b(cell,p,1) = -(vec_a(cell,p,0)*vec_b(cell,p,2)-vec_a(cell,p,2)*vec_b(cell,p,0));
89 vec_a_cross_vec_b(cell,p,2) = vec_a(cell,p,0)*vec_b(cell,p,1)-vec_a(cell,p,1)*vec_b(cell,p,0);
90 }
91 }
92 }
93}
94
95//**********************************************************************
96
97}
98
99#endif
CrossProduct(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
int num_cells
DEPRECATED - use: numCells()