Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Workset_Builder_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_Workset_Builder_impl_hpp__
12#define __Panzer_Workset_Builder_impl_hpp__
13
14#include <iostream>
15#include <vector>
16#include <map>
17#include <algorithm>
18#include "Panzer_Workset.hpp"
19#include "Panzer_CellData.hpp"
20#include "Panzer_BC.hpp"
23
24#include "Phalanx_DataLayout_MDALayout.hpp"
25
26// Intrepid2
27#include "Shards_CellTopology.hpp"
28#include "Intrepid2_DefaultCubatureFactory.hpp"
29#include "Intrepid2_CellTools.hpp"
30#include "Intrepid2_FunctionSpaceTools.hpp"
31#include "Intrepid2_Basis.hpp"
32
33template<typename ArrayT>
34Teuchos::RCP< std::vector<panzer::Workset> >
35panzer::buildWorksets(const WorksetNeeds & needs,
36 const std::string & elementBlock,
37 const std::vector<std::size_t>& local_cell_ids,
38 const ArrayT& node_coordinates)
39{
40 using std::vector;
41 using std::string;
42 using Teuchos::RCP;
43 using Teuchos::rcp;
44
45 panzer::MDFieldArrayFactory mdArrayFactory("",true);
46
47 std::size_t total_num_cells = local_cell_ids.size();
48
49 std::size_t workset_size = needs.cellData.numCells();
50
51 Teuchos::RCP< std::vector<panzer::Workset> > worksets_ptr =
52 Teuchos::rcp(new std::vector<panzer::Workset>);
53 std::vector<panzer::Workset>& worksets = *worksets_ptr;
54
55 // special case for 0 elements!
56 if(total_num_cells==0) {
57
58 // Setup integration rules and basis
59 RCP<vector<int> > ir_degrees = rcp(new vector<int>(0));
60 RCP<vector<string> > basis_names = rcp(new vector<string>(0));
61
62 worksets.resize(1);
63 std::vector<panzer::Workset>::iterator i = worksets.begin();
64 i->setNumberOfCells(0,0,0);
65 i->block_id = elementBlock;
66 i->ir_degrees = ir_degrees;
67 i->basis_names = basis_names;
68
69 for (std::size_t j=0;j<needs.int_rules.size();j++) {
70
71 RCP<panzer::IntegrationValues2<double> > iv2 =
72 rcp(new panzer::IntegrationValues2<double>("",true));
73 iv2->setupArrays(needs.int_rules[j]);
74
75 ir_degrees->push_back(needs.int_rules[j]->cubature_degree);
76 i->int_rules.push_back(iv2);
77 }
78
79 // Need to create all combinations of basis/ir pairings
80 for (std::size_t j=0;j<needs.int_rules.size();j++) {
81 for (std::size_t b=0;b<needs.bases.size();b++) {
82 RCP<panzer::BasisIRLayout> b_layout
83 = rcp(new panzer::BasisIRLayout(needs.bases[b],*needs.int_rules[j]));
84
85 RCP<panzer::BasisValues2<double> > bv2
86 = rcp(new panzer::BasisValues2<double>("",true,true));
87 bv2->setupArrays(b_layout);
88 i->bases.push_back(bv2);
89
90 basis_names->push_back(b_layout->name());
91 }
92
93 }
94
95 return worksets_ptr;
96 } // end special case
97
98 {
99 std::size_t num_worksets = total_num_cells / workset_size;
100 bool last_set_is_full = true;
101 std::size_t last_workset_size = total_num_cells % workset_size;
102 if (last_workset_size != 0) {
103 num_worksets += 1;
104 last_set_is_full = false;
105 }
106
107 worksets.resize(num_worksets);
108 std::vector<panzer::Workset>::iterator i;
109 for (i = worksets.begin(); i != worksets.end(); ++i)
110 i->setNumberOfCells(workset_size,0,0);
111
112 if (!last_set_is_full) {
113 worksets.back().setNumberOfCells(last_workset_size,0,0);
114 }
115 }
116
117 // assign workset cell local ids
118 std::vector<std::size_t>::const_iterator local_begin = local_cell_ids.begin();
119 for (std::vector<panzer::Workset>::iterator wkst = worksets.begin(); wkst != worksets.end(); ++wkst) {
120 std::vector<std::size_t>::const_iterator begin_iter = local_begin;
121 std::vector<std::size_t>::const_iterator end_iter = begin_iter + wkst->num_cells;
122 local_begin = end_iter;
123 wkst->cell_local_ids.assign(begin_iter,end_iter);
124
125 PHX::View<int*> cell_local_ids_k = PHX::View<int*>("Workset:cell_local_ids",wkst->cell_local_ids.size());
126 auto cell_local_ids_k_h = Kokkos::create_mirror_view(cell_local_ids_k);
127 for(std::size_t i=0;i<wkst->cell_local_ids.size();i++)
128 cell_local_ids_k_h(i) = wkst->cell_local_ids[i];
129 Kokkos::deep_copy(cell_local_ids_k, cell_local_ids_k_h);
130 wkst->cell_local_ids_k = cell_local_ids_k;
131
132 wkst->cell_node_coordinates = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("cnc",workset_size,
133 node_coordinates.extent(1),
134 node_coordinates.extent(2));
135 wkst->block_id = elementBlock;
136 wkst->subcell_dim = needs.cellData.baseCellDimension();
137 wkst->subcell_index = 0;
138 }
139
140 TEUCHOS_ASSERT(local_begin == local_cell_ids.end());
141
142 // Copy cell node coordinates into local workset arrays
143 std::size_t offset = 0;
144 for (std::vector<panzer::Workset>::iterator wkst = worksets.begin(); wkst != worksets.end(); ++wkst) {
145 auto cell_node_coordinates = wkst->cell_node_coordinates.get_static_view();
146 Kokkos::parallel_for(wkst->num_cells, KOKKOS_LAMBDA (int cell) {
147 for (std::size_t node = 0; node < node_coordinates.extent(1); ++ node)
148 for (std::size_t dim = 0; dim < node_coordinates.extent(2); ++ dim) {
149 cell_node_coordinates(cell,node,dim) = node_coordinates(cell + offset,node,dim);
150 }
151 });
152 Kokkos::fence();
153 offset += wkst->num_cells;
154 }
155
156 TEUCHOS_ASSERT(offset == Teuchos::as<std::size_t>(node_coordinates.extent(0)));
157
158 // Set ir and basis arrayskset
159 RCP<vector<int> > ir_degrees = rcp(new vector<int>(0));
160 RCP<vector<string> > basis_names = rcp(new vector<string>(0));
161 for (std::vector<panzer::Workset>::iterator wkst = worksets.begin(); wkst != worksets.end(); ++wkst) {
162 wkst->ir_degrees = ir_degrees;
163 wkst->basis_names = basis_names;
164 }
165
166 // setup the integration rules and bases
167 for(std::vector<panzer::Workset>::iterator wkst = worksets.begin(); wkst != worksets.end(); ++wkst)
168 populateValueArrays(wkst->num_cells,false,needs,*wkst);
169
170 return worksets_ptr;
171}
172
173// ****************************************************************
174// ****************************************************************
175
176template<typename ArrayT>
177Teuchos::RCP<std::map<unsigned,panzer::Workset> >
178panzer::buildBCWorkset(const WorksetNeeds & needs,
179 const std::string & elementBlock,
180 const std::vector<std::size_t>& local_cell_ids,
181 const std::vector<std::size_t>& local_side_ids,
182 const ArrayT& node_coordinates,
183 const bool populate_value_arrays)
184{
185 using Teuchos::RCP;
186 using Teuchos::rcp;
187
188 panzer::MDFieldArrayFactory mdArrayFactory("",true);
189
190 // key is local face index, value is workset with all elements
191 // for that local face
192 auto worksets_ptr = Teuchos::rcp(new std::map<unsigned,panzer::Workset>);
193
194 // All elements of boundary condition should go into one workset.
195 // However due to design of Intrepid2 (requires same basis for all
196 // cells), we have to separate the workset based on the local side
197 // index. Each workset for a boundary condition is associated with
198 // a local side for the element
199
200 TEUCHOS_ASSERT(local_side_ids.size() == local_cell_ids.size());
201 TEUCHOS_ASSERT(local_side_ids.size() == static_cast<std::size_t>(node_coordinates.extent(0)));
202
203 // key is local face index, value is a pair of cell index and vector of element local ids
204 std::map< unsigned, std::vector< std::pair< std::size_t, std::size_t>>> element_list;
205 for (std::size_t cell = 0; cell < local_cell_ids.size(); ++cell)
206 element_list[local_side_ids[cell]].push_back(std::make_pair(cell, local_cell_ids[cell]));
207
208 auto& worksets = *worksets_ptr;
209
210 // create worksets
211 for (const auto& side : element_list) {
212
213 auto& cell_local_ids = worksets[side.first].cell_local_ids;
214
215 worksets[side.first].cell_node_coordinates = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("cnc",
216 side.second.size(),
217 node_coordinates.extent(1),
218 node_coordinates.extent(2));
219 auto coords_view = worksets[side.first].cell_node_coordinates.get_view();
220 auto coords_h = Kokkos::create_mirror_view(coords_view);
221
222 auto node_coordinates_h = Kokkos::create_mirror_view(PHX::as_view(node_coordinates));
223 Kokkos::deep_copy(node_coordinates_h, PHX::as_view(node_coordinates));
224
225 for (std::size_t cell = 0; cell < side.second.size(); ++cell) {
226 cell_local_ids.push_back(side.second[cell].second);
227 const auto dim0 = side.second[cell].first;
228
229 for(std::size_t node = 0; node < node_coordinates.extent(1); ++node)
230 {
231 const auto extent = Teuchos::as<std::size_t>(node_coordinates.extent(2));
232
233 for (std::size_t dim = 0; dim < extent; ++dim)
234 coords_h(cell, node, dim) = node_coordinates_h(dim0, node, dim);
235 }
236 }
237
238 Kokkos::deep_copy(coords_view, coords_h);
239
240 const auto cell_local_ids_size = worksets[side.first].cell_local_ids.size();
241 auto cell_local_ids_k = PHX::View<int*>("Workset:cell_local_ids", cell_local_ids_size);
242 auto cell_local_ids_k_h = Kokkos::create_mirror_view(cell_local_ids_k);
243
244 for(std::size_t i = 0; i < cell_local_ids_size; ++i){
245 cell_local_ids_k_h(i) = worksets.at(side.first).cell_local_ids[i];
246 }
247
248 Kokkos::deep_copy(cell_local_ids_k, cell_local_ids_k_h);
249
250 worksets[side.first].cell_local_ids_k = cell_local_ids_k;
251 worksets[side.first].num_cells = worksets[side.first].cell_local_ids.size();
252 worksets[side.first].block_id = elementBlock;
253 worksets[side.first].subcell_dim = needs.cellData.baseCellDimension() - 1;
254 worksets[side.first].subcell_index = side.first;
255 }
256
257 if (populate_value_arrays) {
258 // setup the integration rules and bases
259 for (std::map<unsigned,panzer::Workset>::iterator wkst = worksets.begin();
260 wkst != worksets.end(); ++wkst) {
261
262 populateValueArrays(wkst->second.num_cells,true,needs,wkst->second); // populate "side" values
263 }
264 }
265
266 return worksets_ptr;
267}
268
269// ****************************************************************
270// ****************************************************************
271
272namespace panzer {
273namespace impl {
274
275/* Associate two sets of local side IDs into lists. Each list L has the property
276 * that every local side id in that list is the same, and this holds for each
277 * local side ID set. The smallest set of lists is found. The motivation for
278 * this procedure is to find a 1-1 workset pairing in advance. See the comment
279 * re: Intrepid2 in buildBCWorkset for more.
280 * The return value is an RCP to a map. Only the map's values are of interest
281 * in practice. Each value is a list L. The map's key is a pair (side ID a, side
282 * ID b) that gives rise to the list. We return a pointer to a map so that the
283 * caller does not have to deal with the map type; auto suffices.
284 */
285Teuchos::RCP< std::map<std::pair<std::size_t, std::size_t>, std::vector<std::size_t> > >
286associateCellsBySideIds(const std::vector<std::size_t>& sia /* local_side_ids_a */,
287 const std::vector<std::size_t>& sib /* local_side_ids_b */)
288{
289 TEUCHOS_ASSERT(sia.size() == sib.size());
290
291 auto sip2i_p = Teuchos::rcp(new std::map< std::pair<std::size_t, std::size_t>, std::vector<std::size_t> >);
292 auto& sip2i = *sip2i_p;
293
294 for (std::size_t i = 0; i < sia.size(); ++i)
295 sip2i[std::make_pair(sia[i], sib[i])].push_back(i);
296
297 return sip2i_p;
298}
299
300// Set s = a(idxs). No error checking.
301template <typename T>
302void subset(const std::vector<T>& a, const std::vector<std::size_t>& idxs, std::vector<T>& s)
303{
304 s.resize(idxs.size());
305 for (std::size_t i = 0; i < idxs.size(); ++i)
306 s[i] = a[idxs[i]];
307}
308
309template<typename ArrayT>
310Teuchos::RCP<std::map<unsigned,panzer::Workset> >
312 const std::string & blockid_a,
313 const std::vector<std::size_t>& local_cell_ids_a,
314 const std::vector<std::size_t>& local_side_ids_a,
315 const ArrayT& node_coordinates_a,
316 const panzer::WorksetNeeds & needs_b,
317 const std::string & blockid_b,
318 const std::vector<std::size_t>& local_cell_ids_b,
319 const std::vector<std::size_t>& local_side_ids_b,
320 const ArrayT& node_coordinates_b,
321 const WorksetNeeds& needs_b2)
322{
323 TEUCHOS_ASSERT(local_cell_ids_a.size() == local_cell_ids_b.size());
324 // Get a and b workset maps separately, but don't populate b's arrays.
325 const Teuchos::RCP<std::map<unsigned,panzer::Workset> >
326 mwa = buildBCWorkset(needs_a,blockid_a, local_cell_ids_a, local_side_ids_a, node_coordinates_a),
327 mwb = buildBCWorkset(needs_b2, blockid_b, local_cell_ids_b, local_side_ids_b,
328 node_coordinates_b, false /* populate_value_arrays */);
329 TEUCHOS_ASSERT(mwa->size() == 1 && mwb->size() == 1);
330 for (std::map<unsigned,panzer::Workset>::iterator ait = mwa->begin(), bit = mwb->begin();
331 ait != mwa->end(); ++ait, ++bit) {
332 TEUCHOS_ASSERT(Teuchos::as<std::size_t>(ait->second.num_cells) == local_cell_ids_a.size() &&
333 Teuchos::as<std::size_t>(bit->second.num_cells) == local_cell_ids_b.size());
334 panzer::Workset& wa = ait->second;
335 // Copy b's details(0) to a's details(1).
336 wa.other = Teuchos::rcp(new panzer::WorksetDetails(bit->second.details(0)));
337 // Populate details(1) arrays so that IP are in order corresponding to details(0).
338 populateValueArrays(wa.num_cells, true, needs_b, wa.details(1), Teuchos::rcpFromRef(wa.details(0)));
339 }
340 // Now mwa has everything we need.
341 return mwa;
342}
343
344} // namespace impl
345} // namespace panzer
346
347// ****************************************************************
348// ****************************************************************
349
350template<typename ArrayT>
351Teuchos::RCP<std::map<unsigned,panzer::Workset> >
353 const std::string & blockid_a,
354 const std::vector<std::size_t>& local_cell_ids_a,
355 const std::vector<std::size_t>& local_side_ids_a,
356 const ArrayT& node_coordinates_a,
357 const panzer::WorksetNeeds & needs_b,
358 const std::string & blockid_b,
359 const std::vector<std::size_t>& local_cell_ids_b,
360 const std::vector<std::size_t>& local_side_ids_b,
361 const ArrayT& node_coordinates_b)
362{
363 // Since Intrepid2 requires all side IDs in a workset to be the same (see
364 // Intrepid2 comment above), we break the element list into pieces such that
365 // each piece contains elements on each side of the interface L_a and L_b and
366 // all elemnets L_a have the same side ID, and the same for L_b.
367 auto side_id_associations = impl::associateCellsBySideIds(local_side_ids_a, local_side_ids_b);
368 if (side_id_associations->size() == 1) {
369 // Common case of one workset on each side; optimize for it.
370 return impl::buildBCWorksetForUniqueSideId(needs_a, blockid_a, local_cell_ids_a, local_side_ids_a, node_coordinates_a,
371 needs_b, blockid_b, local_cell_ids_b, local_side_ids_b, node_coordinates_b,
372 needs_b);
373 } else {
374 // The interface has elements having a mix of side IDs, so deal with each
375 // pair in turn.
376 Teuchos::RCP<std::map<unsigned,panzer::Workset> > mwa = Teuchos::rcp(new std::map<unsigned,panzer::Workset>);
377 std::vector<std::size_t> lci_a, lci_b, lsi_a, lsi_b;
378 panzer::MDFieldArrayFactory mdArrayFactory("", true);
379 const int d1 = Teuchos::as<std::size_t>(node_coordinates_a.extent(1)),
380 d2 = Teuchos::as<std::size_t>(node_coordinates_a.extent(2));
381 for (auto it : *side_id_associations) {
382 const auto& idxs = it.second;
383 impl::subset(local_cell_ids_a, idxs, lci_a);
384 impl::subset(local_side_ids_a, idxs, lsi_a);
385 impl::subset(local_cell_ids_b, idxs, lci_b);
386 impl::subset(local_side_ids_b, idxs, lsi_b);
387 auto nc_a = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("nc_a", idxs.size(), d1, d2);
388 auto nc_b = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("nc_b", idxs.size(), d1, d2);
389 auto nc_a_h = Kokkos::create_mirror_view(nc_a.get_static_view());
390 auto nc_b_h = Kokkos::create_mirror_view(nc_b.get_static_view());
391 auto node_coordinates_a_h = Kokkos::create_mirror_view(PHX::as_view(node_coordinates_a));
392 auto node_coordinates_b_h = Kokkos::create_mirror_view(PHX::as_view(node_coordinates_b));
393 Kokkos::deep_copy(node_coordinates_a_h, PHX::as_view(node_coordinates_a));
394 Kokkos::deep_copy(node_coordinates_b_h, PHX::as_view(node_coordinates_b));
395 for (std::size_t i = 0; i < idxs.size(); ++i) {
396 const auto ii = idxs[i];
397 for (int j = 0; j < d1; ++j)
398 for (int k = 0; k < d2; ++k) {
399 nc_a_h(i, j, k) = node_coordinates_a_h(ii, j, k);
400 nc_b_h(i, j, k) = node_coordinates_b_h(ii, j, k);
401 }
402 }
403 Kokkos::deep_copy(nc_a.get_static_view(), nc_a_h);
404 Kokkos::deep_copy(nc_b.get_static_view(), nc_b_h);
405 auto mwa_it = impl::buildBCWorksetForUniqueSideId(needs_a,blockid_a, lci_a, lsi_a, nc_a,
406 needs_b,blockid_b, lci_b, lsi_b, nc_b,
407 needs_b);
408 TEUCHOS_ASSERT(mwa_it->size() == 1);
409 // Form a unique key that encodes the pair (side ID a, side ID b). We
410 // abuse the key here in the sense that it is everywhere else understood
411 // to correspond to the side ID of the elements in the workset.
412 // 1000 is a number substantially larger than is needed for any element.
413 const std::size_t key = lsi_a[0] * 1000 + lsi_b[0];
414 (*mwa)[key] = mwa_it->begin()->second;
415 }
416 return mwa;
417 }
418}
419
420// ****************************************************************
421// ****************************************************************
422
423template<typename ArrayT>
424Teuchos::RCP<std::vector<panzer::Workset> >
425panzer::buildEdgeWorksets(const WorksetNeeds & needs_a,
426 const std::string & eblock_a,
427 const std::vector<std::size_t>& local_cell_ids_a,
428 const std::vector<std::size_t>& local_side_ids_a,
429 const ArrayT& node_coordinates_a,
430 const WorksetNeeds & needs_b,
431 const std::string & eblock_b,
432 const std::vector<std::size_t>& local_cell_ids_b,
433 const std::vector<std::size_t>& local_side_ids_b,
434 const ArrayT& node_coordinates_b)
435{
436 using Teuchos::RCP;
437 using Teuchos::rcp;
438
439 panzer::MDFieldArrayFactory mdArrayFactory("",true);
440
441 std::size_t total_num_cells_a = local_cell_ids_a.size();
442 std::size_t total_num_cells_b = local_cell_ids_b.size();
443
444 TEUCHOS_ASSERT(total_num_cells_a==total_num_cells_b);
445 TEUCHOS_ASSERT(local_side_ids_a.size() == local_cell_ids_a.size());
446 TEUCHOS_ASSERT(local_side_ids_a.size() == static_cast<std::size_t>(node_coordinates_a.extent(0)));
447 TEUCHOS_ASSERT(local_side_ids_b.size() == local_cell_ids_b.size());
448 TEUCHOS_ASSERT(local_side_ids_b.size() == static_cast<std::size_t>(node_coordinates_b.extent(0)));
449
450 std::size_t total_num_cells = total_num_cells_a;
451
452 std::size_t workset_size = needs_a.cellData.numCells();
453
454 Teuchos::RCP< std::vector<panzer::Workset> > worksets_ptr =
455 Teuchos::rcp(new std::vector<panzer::Workset>);
456 std::vector<panzer::Workset>& worksets = *worksets_ptr;
457
458 // special case for 0 elements!
459 if(total_num_cells==0) {
460
461 // Setup integration rules and basis
462 RCP<std::vector<int> > ir_degrees = rcp(new std::vector<int>(0));
463 RCP<std::vector<std::string> > basis_names = rcp(new std::vector<std::string>(0));
464
465 worksets.resize(1);
466 std::vector<panzer::Workset>::iterator i = worksets.begin();
467
468 i->details(0).block_id = eblock_a;
469 i->other = Teuchos::rcp(new panzer::WorksetDetails);
470 i->details(1).block_id = eblock_b;
471
472 i->num_cells = 0;
473 i->ir_degrees = ir_degrees;
474 i->basis_names = basis_names;
475
476 for(std::size_t j=0;j<needs_a.int_rules.size();j++) {
477
478 RCP<panzer::IntegrationValues2<double> > iv2 =
479 rcp(new panzer::IntegrationValues2<double>("",true));
480 iv2->setupArrays(needs_a.int_rules[j]);
481
482 ir_degrees->push_back(needs_a.int_rules[j]->cubature_degree);
483 i->int_rules.push_back(iv2);
484 }
485
486 // Need to create all combinations of basis/ir pairings
487 for(std::size_t j=0;j<needs_a.int_rules.size();j++) {
488
489 for(std::size_t b=0;b<needs_a.bases.size();b++) {
490
491 RCP<panzer::BasisIRLayout> b_layout = rcp(new panzer::BasisIRLayout(needs_a.bases[b],*needs_a.int_rules[j]));
492
493 RCP<panzer::BasisValues2<double> > bv2 =
494 rcp(new panzer::BasisValues2<double>("",true,true));
495 bv2->setupArrays(b_layout);
496 i->bases.push_back(bv2);
497
498 basis_names->push_back(b_layout->name());
499 }
500
501 }
502
503 return worksets_ptr;
504 } // end special case
505
506 // This collects all the elements that share the same sub cell pairs, this makes it easier to
507 // build the required worksets
508 // key is the pair of local face indices, value is a vector of cell indices that satisfy this pair
509 std::map<std::pair<unsigned,unsigned>,std::vector<std::size_t> > element_list;
510 for (std::size_t cell=0; cell < local_cell_ids_a.size(); ++cell)
511 element_list[std::make_pair(local_side_ids_a[cell],local_side_ids_b[cell])].push_back(cell);
512
513 // figure out how many worksets will be needed, resize workset vector accordingly
514 std::size_t num_worksets = 0;
515 for(const auto& edge : element_list) {
516 std::size_t num_worksets_for_edge = edge.second.size() / workset_size;
517 std::size_t last_workset_size = edge.second.size() % workset_size;
518 if(last_workset_size!=0)
519 num_worksets_for_edge += 1;
520
521 num_worksets += num_worksets_for_edge;
522 }
523 worksets.resize(num_worksets);
524
525 // fill the worksets
526 std::vector<Workset>::iterator current_workset = worksets.begin();
527 for(const auto& edge : element_list) {
528 // loop over each workset
529 const std::vector<std::size_t> & cell_indices = edge.second;
530
531 current_workset = buildEdgeWorksets(cell_indices,
532 needs_a,eblock_a,local_cell_ids_a,local_side_ids_a,node_coordinates_a,
533 needs_b,eblock_b,local_cell_ids_b,local_side_ids_b,node_coordinates_b,
534 current_workset);
535 }
536
537 // sanity check
538 TEUCHOS_ASSERT(current_workset==worksets.end());
539
540 return worksets_ptr;
541}
542
543template<typename ArrayT>
544std::vector<panzer::Workset>::iterator
545panzer::buildEdgeWorksets(const std::vector<std::size_t> & cell_indices,
546 const WorksetNeeds & needs_a,
547 const std::string & eblock_a,
548 const std::vector<std::size_t>& local_cell_ids_a,
549 const std::vector<std::size_t>& local_side_ids_a,
550 const ArrayT& node_coordinates_a,
551 const WorksetNeeds & needs_b,
552 const std::string & eblock_b,
553 const std::vector<std::size_t>& local_cell_ids_b,
554 const std::vector<std::size_t>& local_side_ids_b,
555 const ArrayT& node_coordinates_b,
556 std::vector<Workset>::iterator beg)
557{
558 panzer::MDFieldArrayFactory mdArrayFactory("",true);
559
560 std::vector<Workset>::iterator wkst = beg;
561
562 std::size_t current_cell_index = 0;
563 while (current_cell_index<cell_indices.size()) {
564 std::size_t workset_size = needs_a.cellData.numCells();
565
566 // allocate workset details (details(0) is already associated with the
567 // workset object itself)
568 wkst->other = Teuchos::rcp(new panzer::WorksetDetails);
569
570 wkst->subcell_dim = needs_a.cellData.baseCellDimension()-1;
571
572 wkst->details(0).subcell_index = local_side_ids_a[cell_indices[current_cell_index]];
573 wkst->details(0).block_id = eblock_a;
574 wkst->details(0).cell_node_coordinates = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("cnc",workset_size,
575 node_coordinates_a.extent(1),
576 node_coordinates_a.extent(2));
577
578 wkst->details(1).subcell_index = local_side_ids_b[cell_indices[current_cell_index]];
579 wkst->details(1).block_id = eblock_b;
580 wkst->details(1).cell_node_coordinates = mdArrayFactory.buildStaticArray<double,Cell,NODE,Dim>("cnc",workset_size,
581 node_coordinates_a.extent(1),
582 node_coordinates_a.extent(2));
583
584 std::size_t remaining_cells = cell_indices.size()-current_cell_index;
585 if(remaining_cells<workset_size)
586 workset_size = remaining_cells;
587
588 // this is the true number of cells in this workset
589 wkst->setNumberOfCells(workset_size,0,0);
590 wkst->details(0).cell_local_ids.resize(workset_size);
591 wkst->details(1).cell_local_ids.resize(workset_size);
592
593 auto dim0_cell_node_coordinates_view = wkst->details(0).cell_node_coordinates.get_static_view();
594 auto dim0_cell_node_coordinates_h = Kokkos::create_mirror_view(dim0_cell_node_coordinates_view);
595 Kokkos::deep_copy(dim0_cell_node_coordinates_h, dim0_cell_node_coordinates_view);
596
597 auto dim1_cell_node_coordinates_view = wkst->details(1).cell_node_coordinates.get_static_view();
598 auto dim1_cell_node_coordinates_h = Kokkos::create_mirror_view(dim1_cell_node_coordinates_view);
599 Kokkos::deep_copy(dim1_cell_node_coordinates_h, dim1_cell_node_coordinates_view);
600
601 auto node_coordinates_a_h = Kokkos::create_mirror_view(node_coordinates_a);
602 Kokkos::deep_copy(node_coordinates_a_h, node_coordinates_a);
603
604 auto node_coordinates_b_h = Kokkos::create_mirror_view(node_coordinates_b);
605 Kokkos::deep_copy(node_coordinates_b_h, node_coordinates_b);
606
607 for(std::size_t cell=0;cell<workset_size; cell++,current_cell_index++) {
608
609 wkst->details(0).cell_local_ids[cell] = local_cell_ids_a[cell_indices[current_cell_index]];
610 wkst->details(1).cell_local_ids[cell] = local_cell_ids_b[cell_indices[current_cell_index]];
611
612 for (std::size_t node = 0; node < Teuchos::as<std::size_t>(node_coordinates_a.extent(1)); ++ node) {
613 for (std::size_t dim = 0; dim < Teuchos::as<std::size_t>(node_coordinates_a.extent(2)); ++ dim) {
614 dim0_cell_node_coordinates_h(cell,node,dim) = node_coordinates_a_h(cell_indices[current_cell_index],node,dim);
615 dim1_cell_node_coordinates_h(cell,node,dim) = node_coordinates_b_h(cell_indices[current_cell_index],node,dim);
616 }
617 }
618 }
619
620 Kokkos::deep_copy(dim0_cell_node_coordinates_view, dim0_cell_node_coordinates_h);
621 Kokkos::deep_copy(dim1_cell_node_coordinates_view, dim1_cell_node_coordinates_h);
622
623 auto cell_local_ids_k_0 = PHX::View<int*>("Workset:cell_local_ids",wkst->details(0).cell_local_ids.size());
624 auto cell_local_ids_k_0_h = Kokkos::create_mirror_view(cell_local_ids_k_0);
625
626 auto cell_local_ids_k_1 = PHX::View<int*>("Workset:cell_local_ids",wkst->details(1).cell_local_ids.size());
627 auto cell_local_ids_k_1_h = Kokkos::create_mirror_view(cell_local_ids_k_1);
628
629 for(std::size_t i=0;i<wkst->details(0).cell_local_ids.size();i++)
630 cell_local_ids_k_0_h(i) = wkst->details(0).cell_local_ids[i];
631 for(std::size_t i=0;i<wkst->details(1).cell_local_ids.size();i++)
632 cell_local_ids_k_1_h(i) = wkst->details(1).cell_local_ids[i];
633
634 Kokkos::deep_copy(cell_local_ids_k_0, cell_local_ids_k_0_h);
635 Kokkos::deep_copy(cell_local_ids_k_1, cell_local_ids_k_1_h);
636
637 wkst->details(0).cell_local_ids_k = cell_local_ids_k_0;
638 wkst->details(1).cell_local_ids_k = cell_local_ids_k_1;
639
640 // fill the BasisValues and IntegrationValues arrays
641 std::size_t max_workset_size = needs_a.cellData.numCells();
642 populateValueArrays(max_workset_size,true,needs_a,wkst->details(0)); // populate "side" values
643 populateValueArrays(max_workset_size,true,needs_b,wkst->details(1),Teuchos::rcpFromRef(wkst->details(0)));
644
645 wkst++;
646 }
647
648 return wkst;
649}
650
651#endif
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
int num_cells
DEPRECATED - use: numCells()
WorksetDetails & details(const int i)
Convenience wrapper to operator() for pointer access.
Teuchos::RCP< WorksetDetails > other
Teuchos::RCP< std::map< unsigned, panzer::Workset > > buildBCWorksetForUniqueSideId(const panzer::WorksetNeeds &needs_a, const std::string &blockid_a, const std::vector< std::size_t > &local_cell_ids_a, const std::vector< std::size_t > &local_side_ids_a, const ArrayT &node_coordinates_a, const panzer::WorksetNeeds &needs_b, const std::string &blockid_b, const std::vector< std::size_t > &local_cell_ids_b, const std::vector< std::size_t > &local_side_ids_b, const ArrayT &node_coordinates_b, const WorksetNeeds &needs_b2)
void subset(const std::vector< T > &a, const std::vector< std::size_t > &idxs, std::vector< T > &s)
Teuchos::RCP< std::map< std::pair< std::size_t, std::size_t >, std::vector< std::size_t > > > associateCellsBySideIds(const std::vector< std::size_t > &sia, const std::vector< std::size_t > &sib)
void populateValueArrays(std::size_t num_cells, bool isSide, const WorksetNeeds &needs, WorksetDetails &details, const Teuchos::RCP< WorksetDetails > other_details)
Teuchos::RCP< std::vector< Workset > > buildWorksets(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const ArrayT &node_coordinates)
Teuchos::RCP< std::vector< Workset > > buildEdgeWorksets(const WorksetNeeds &needs_a, const std::string &eblock_a, const std::vector< std::size_t > &local_cell_ids_a, const std::vector< std::size_t > &local_side_ids_a, const ArrayT &node_coordinates_a, const WorksetNeeds &needs_b, const std::string &eblock_b, const std::vector< std::size_t > &local_cell_ids_b, const std::vector< std::size_t > &local_side_ids_b, const ArrayT &node_coordinates_b)
Teuchos::RCP< std::map< unsigned, Workset > > buildBCWorkset(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const std::vector< std::size_t > &local_side_ids, const ArrayT &node_coordinates, const bool populate_value_arrays=true)