Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_FaceToElement.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/*
12 * FaceToElement.hpp
13 *
14 * Created on: Nov 15, 2016
15 * Author: mbetten
16 */
17
18#ifndef PANZER_FACE_TO_ELEMENT_HPP
19#define PANZER_FACE_TO_ELEMENT_HPP
20
21#include "Phalanx_KokkosDeviceTypes.hpp"
22
24
25#include <Tpetra_Map.hpp>
26#include <Tpetra_MultiVector.hpp>
27#include <Tpetra_Import.hpp>
28#include <Tpetra_Export.hpp>
29
30namespace panzer
31{
32
36template <typename LocalOrdinal,typename GlobalOrdinal>
38private:
39 FaceToElement(const FaceToElement &); // disallowed
40
41public:
42
44
46 const Teuchos::RCP<const Teuchos::Comm<int>> comm);
47
51 const Teuchos::RCP<const Teuchos::Comm<int>> comm);
52
53
54 GlobalOrdinal getLeftElem (GlobalOrdinal face_id) const
55 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return elems_by_face_(lid,0);}
56
57 GlobalOrdinal getRightElem(GlobalOrdinal face_id) const
58 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return elems_by_face_(lid,1);}
59
60 int getLeftBlock (GlobalOrdinal face_id) const
61 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return blocks_by_face_(lid,0);}
62
63 int getRightBlock(GlobalOrdinal face_id) const
64 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return blocks_by_face_(lid,1);}
65
66 int getLeftProc (GlobalOrdinal face_id) const
67 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return procs_by_face_(lid,0);}
68
69 int getRightProc (GlobalOrdinal face_id) const
70 {LocalOrdinal lid = face_map_->getLocalElement(face_id); return procs_by_face_(lid,1);}
71
72 PHX::View<const GlobalOrdinal*[2]> getFaceToElementsMap() const
73 { return elems_by_face_; }
74
75 PHX::View<const int*[2]> getFaceToCellLocalIdxMap() const
76 { return lidx_by_face_; }
77
78protected:
79
80 PHX::View<GlobalOrdinal *[2]> elems_by_face_;
81 PHX::View<int *[2]> lidx_by_face_;
82 PHX::View<int *[2]> blocks_by_face_;
83 PHX::View<int *[2]> procs_by_face_;
84
85 typedef Tpetra::KokkosCompat::KokkosDeviceWrapperNode<PHX::Device> NodeType;
86 typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, NodeType> Map;
87 typedef Tpetra::Export<LocalOrdinal, GlobalOrdinal, NodeType> Export;
88 typedef Tpetra::Import<LocalOrdinal, GlobalOrdinal, NodeType> Import;
89 typedef Tpetra::MultiVector<GlobalOrdinal, LocalOrdinal, GlobalOrdinal, NodeType> GOMultiVector;
90
91
92 Teuchos::RCP<const Map> face_map_;
93
94};
95
96}
97
98#endif
Pure virtual base class for supplying mesh connectivity information to the DOF Manager.
Teuchos::RCP< const Map > face_map_
Tpetra::MultiVector< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, NodeType > GOMultiVector
Tpetra::Export< LocalOrdinal, GlobalOrdinal, NodeType > Export
PHX::View< int *[2]> lidx_by_face_
GlobalOrdinal getLeftElem(GlobalOrdinal face_id) const
void initialize(panzer::ConnManager &conn, const Teuchos::RCP< const Teuchos::Comm< int > > comm)
PHX::View< GlobalOrdinal *[2]> elems_by_face_
GlobalOrdinal getRightElem(GlobalOrdinal face_id) const
int getRightProc(GlobalOrdinal face_id) const
int getRightBlock(GlobalOrdinal face_id) const
int getLeftProc(GlobalOrdinal face_id) const
PHX::View< int *[2]> procs_by_face_
PHX::View< const int *[2]> getFaceToCellLocalIdxMap() const
Tpetra::Map< LocalOrdinal, GlobalOrdinal, NodeType > Map
Tpetra::Import< LocalOrdinal, GlobalOrdinal, NodeType > Import
Tpetra::KokkosCompat::KokkosDeviceWrapperNode< PHX::Device > NodeType
int getLeftBlock(GlobalOrdinal face_id) const
PHX::View< int *[2]> blocks_by_face_
PHX::View< const GlobalOrdinal *[2]> getFaceToElementsMap() const
FaceToElement(const FaceToElement &)