18 #ifndef PANZER_FACE_TO_ELEMENT_IMPL_HPP
19 #define PANZER_FACE_TO_ELEMENT_IMPL_HPP
34 template <
typename LocalOrdinal,
typename GlobalOrdinal>
40 #ifndef PANZER_HIDE_DEPRECATED_CODE
45 template <
typename LocalOrdinal,
typename GlobalOrdinal>
53 template <
typename LocalOrdinal,
typename GlobalOrdinal>
58 initialize(conn, comm);
61 #ifndef PANZER_HIDE_DEPRECATED_CODE
66 template <
typename LocalOrdinal,
typename GlobalOrdinal>
72 initialize(conn, comm_world);
76 template <
typename LocalOrdinal,
typename GlobalOrdinal>
83 std::vector<std::string> block_ids;
89 std::vector<shards::CellTopology> ebt;
91 dimension = ebt[0].getDimension();
93 int my_rank = comm->getRank();
95 int nprocs = comm->getSize();
98 std::vector<GlobalOrdinal> element_GIDS;
99 for (
size_t iblk = 0 ; iblk < block_ids.size(); ++iblk) {
102 if ( dimension == 1 ) {
105 }
else if ( dimension == 2 ){
113 const std::vector<LocalOrdinal> &block_elems = conn.
getElementBlock(block_ids[iblk]);
114 for (
size_t i=0; i<block_elems.size(); ++i) {
116 element_GIDS.push_back(*connectivity);
125 std::vector<GlobalOrdinal> face_GIDS;
126 std::set<GlobalOrdinal> set_of_face_GIDS;
127 for (
size_t iblk = 0 ; iblk < block_ids.size(); ++iblk) {
130 if ( dimension == 1 ) {
133 }
else if ( dimension == 2 ){
141 const std::vector<LocalOrdinal> &block_elems = conn.
getElementBlock(block_ids[iblk]);
142 for (
size_t i=0; i<block_elems.size(); ++i) {
144 const panzer::GlobalOrdinal * connectivity = conn.
getConnectivity(block_elems[i]);
145 for (
int iface=0; iface<n_conn; ++iface)
146 set_of_face_GIDS.insert(connectivity[iface]);
149 face_GIDS.insert(face_GIDS.begin(), set_of_face_GIDS.begin(), set_of_face_GIDS.end());
152 owned_face_map = Tpetra::createOneToOne(face_map);
161 face2elem_mv->putScalar(-1-shift);
163 auto f2e = face2elem_mv->getLocalViewHost(Tpetra::Access::ReadWrite);
164 auto b1 = Kokkos::subview(f2e,Kokkos::ALL(),0);
165 auto e1 = Kokkos::subview(f2e,Kokkos::ALL(),1);
166 auto p1 = Kokkos::subview(f2e,Kokkos::ALL(),2);
167 auto l1 = Kokkos::subview(f2e,Kokkos::ALL(),3);
168 auto b2 = Kokkos::subview(f2e,Kokkos::ALL(),4);
169 auto e2 = Kokkos::subview(f2e,Kokkos::ALL(),5);
170 auto p2 = Kokkos::subview(f2e,Kokkos::ALL(),6);
171 auto l2 = Kokkos::subview(f2e,Kokkos::ALL(),7);
173 GlobalOrdinal my_elem = 0;
174 for (
size_t iblk = 0 ; iblk < block_ids.size(); ++iblk) {
177 if ( dimension == 1 ) {
180 }
else if ( dimension == 2 ){
188 const std::vector<LocalOrdinal> &block_elems = conn.
getElementBlock(block_ids[iblk]);
189 for (
size_t i=0; i<block_elems.size(); ++i) {
191 const panzer::GlobalOrdinal * connectivity = conn.
getConnectivity(block_elems[i]);
192 for (
int iface=0; iface<n_conn; ++iface) {
193 LocalOrdinal f = face_map->getLocalElement(connectivity[iface]);
199 e1[f] = elem_map->getGlobalElement(my_elem)+shift;
200 p1[f] = my_rank+shift;
202 }
else if (b2[f] < 0){
204 e2[f] = elem_map->getGlobalElement(my_elem)+shift;
205 p2[f] = my_rank+shift;
216 Import imp(owned_face_map, face_map);
217 Export exp(face_map, owned_face_map);
218 owned_face2elem_mv->doExport(*face2elem_mv, exp, Tpetra::ADD);
221 auto f2e = face2elem_mv->getLocalViewHost(Tpetra::Access::ReadWrite);
222 auto b1 = Kokkos::subview(f2e,Kokkos::ALL(),0);
223 auto e1 = Kokkos::subview(f2e,Kokkos::ALL(),1);
224 auto p1 = Kokkos::subview(f2e,Kokkos::ALL(),2);
225 auto l1 = Kokkos::subview(f2e,Kokkos::ALL(),3);
226 auto b2 = Kokkos::subview(f2e,Kokkos::ALL(),4);
227 auto e2 = Kokkos::subview(f2e,Kokkos::ALL(),5);
228 auto p2 = Kokkos::subview(f2e,Kokkos::ALL(),6);
229 auto l2 = Kokkos::subview(f2e,Kokkos::ALL(),7);
231 auto of2e = owned_face2elem_mv->getLocalViewHost(Tpetra::Access::ReadWrite);
232 auto ob1 = Kokkos::subview(of2e,Kokkos::ALL(),0);
233 auto oe1 = Kokkos::subview(of2e,Kokkos::ALL(),1);
234 auto op1 = Kokkos::subview(of2e,Kokkos::ALL(),2);
235 auto ol1 = Kokkos::subview(of2e,Kokkos::ALL(),3);
236 auto ob2 = Kokkos::subview(of2e,Kokkos::ALL(),4);
237 auto oe2 = Kokkos::subview(of2e,Kokkos::ALL(),5);
238 auto op2 = Kokkos::subview(of2e,Kokkos::ALL(),6);
239 auto ol2 = Kokkos::subview(of2e,Kokkos::ALL(),7);
244 for (
size_t i=0; i<ob1.size();++i){
247 assert(b1[i] >= shift);
249 LocalOrdinal shared_local_id = face_map->getLocalElement(owned_face_map->getGlobalElement(i));
251 if (ob1[i] == b1[shared_local_id] && ob2[i] == b2[shared_local_id] &&
252 oe1[i] == e1[shared_local_id] && oe2[i] == e2[shared_local_id]) {
259 if (ob1[i] < b1[shared_local_id] || oe1[i] < e1[shared_local_id]) {
263 if ( ob1[i] > b1[shared_local_id] || oe1[i] > e1[shared_local_id]) {
265 assert(ob2[i] < 0 && oe2[i] < 0);
266 ob2[i] = ob1[i] - b1[shared_local_id];
267 oe2[i] = oe1[i] - e1[shared_local_id];
268 op2[i] = op1[i] - p1[shared_local_id];
269 ol2[i] = ol1[i] - l1[shared_local_id];
271 ob1[i] = b1[shared_local_id];
272 oe1[i] = e1[shared_local_id];
273 op1[i] = p1[shared_local_id];
274 ol1[i] = l1[shared_local_id];
276 assert(op1[i] >=0 && op2[i] >= 0 && op1[i] < nprocs+shift && op2[i] < nprocs+shift);
280 face2elem_mv->doImport(*owned_face2elem_mv, imp, Tpetra::REPLACE);
285 face_map_ = face_map;
286 LocalOrdinal nfaces = face_map_->getLocalNumElements();
287 elems_by_face_ = PHX::View<GlobalOrdinal *[2]>(
"FaceToElement::elems_by_face_", nfaces);
288 lidx_by_face_ = PHX::View<int *[2]>(
"FaceToElement::elems_by_face_", nfaces);
289 blocks_by_face_ = PHX::View<int *[2]> (
"FaceToElement::blocks_by_face_", nfaces);
290 procs_by_face_ = PHX::View<int *[2]> (
"FaceToElement::procs_by_face_", nfaces);
291 auto elems_by_face_h = Kokkos::create_mirror_view(elems_by_face_);
292 auto blocks_by_face_h = Kokkos::create_mirror_view(blocks_by_face_);
293 auto procs_by_face_h = Kokkos::create_mirror_view(procs_by_face_);
294 auto lidx_by_face_h = Kokkos::create_mirror_view(lidx_by_face_);
296 auto f2e = face2elem_mv->getLocalViewHost(Tpetra::Access::ReadWrite);
297 auto b1 = Kokkos::subview(f2e,Kokkos::ALL(),0);
298 auto e1 = Kokkos::subview(f2e,Kokkos::ALL(),1);
299 auto p1 = Kokkos::subview(f2e,Kokkos::ALL(),2);
300 auto l1 = Kokkos::subview(f2e,Kokkos::ALL(),3);
301 auto b2 = Kokkos::subview(f2e,Kokkos::ALL(),4);
302 auto e2 = Kokkos::subview(f2e,Kokkos::ALL(),5);
303 auto p2 = Kokkos::subview(f2e,Kokkos::ALL(),6);
304 auto l2 = Kokkos::subview(f2e,Kokkos::ALL(),7);
306 for (LocalOrdinal i=0; i< nfaces; ++i) {
307 elems_by_face_h (i,0) = e1[i]-shift;
308 elems_by_face_h (i,1) = e2[i]-shift;
309 blocks_by_face_h(i,0) = b1[i]-shift;
310 blocks_by_face_h(i,1) = b2[i]-shift;
311 procs_by_face_h (i,0) = p1[i]-shift;
312 procs_by_face_h (i,1) = p2[i]-shift;
313 lidx_by_face_h (i,0) = l1[i]-shift;
314 lidx_by_face_h (i,1) = l2[i]-shift;
316 if(elems_by_face_h (i,0) < 0){
317 elems_by_face_h (i,0) = -1;
319 if(elems_by_face_h (i,1) < 0){
320 elems_by_face_h (i,1) = -1;
323 Kokkos::deep_copy(elems_by_face_, elems_by_face_h);
324 Kokkos::deep_copy(blocks_by_face_, blocks_by_face_h);
325 Kokkos::deep_copy(procs_by_face_, procs_by_face_h);
326 Kokkos::deep_copy(lidx_by_face_, lidx_by_face_h);
Tpetra::Import< LocalOrdinal, GlobalOrdinal, NodeType > Import
virtual LocalOrdinal getConnectivitySize(LocalOrdinal localElmtId) const =0
void initialize(panzer::ConnManager &conn)
Tpetra::Map< LocalOrdinal, GlobalOrdinal, NodeType > Map
virtual const std::vector< LocalOrdinal > & getElementBlock(const std::string &blockID) const =0
virtual void getElementBlockIds(std::vector< std::string > &elementBlockIds) const =0
Pure virtual base class for supplying mesh connectivity information to the DOF Manager.
virtual void getElementBlockTopologies(std::vector< shards::CellTopology > &elementBlockTopologies) const =0
Tpetra::Export< LocalOrdinal, GlobalOrdinal, NodeType > Export
Tpetra::MultiVector< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, NodeType > GOMultiVector
virtual void buildConnectivity(const FieldPattern &fp)=0
virtual const GlobalOrdinal * getConnectivity(LocalOrdinal localElmtId) const =0