MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_IntrepidPCoarsenFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_IPCFACTORY_DEF_HPP
47#define MUELU_IPCFACTORY_DEF_HPP
48
49#include <Xpetra_Matrix.hpp>
50#include <Xpetra_IO.hpp>
51#include <sstream>
52#include <algorithm>
53
55
56#include "MueLu_Level.hpp"
57#include "MueLu_MasterList.hpp"
58#include "MueLu_Monitor.hpp"
59#include "MueLu_PerfUtils.hpp"
60#include "MueLu_Utilities.hpp"
61
62#include "Teuchos_ScalarTraits.hpp"
63
64// Intrepid Headers
65
66//Intrepid_HGRAD_HEX_C1_FEM.hpp
67//Intrepid_HGRAD_HEX_C2_FEM.hpp
68//Intrepid_HGRAD_HEX_Cn_FEM.hpp
69//Intrepid_HGRAD_HEX_I2_FEM.hpp
70#include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
71#include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp"
72//Intrepid_HGRAD_LINE_Cn_FEM_JACOBI.hpp
73//Intrepid_HGRAD_POLY_C1_FEM.hpp
74//Intrepid_HGRAD_PYR_C1_FEM.hpp
75//Intrepid_HGRAD_PYR_I2_FEM.hpp
76#include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
77//#include Intrepid_HGRAD_QUAD_C2_FEM.hpp
78#include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp"
79//Intrepid_HGRAD_TET_C1_FEM.hpp
80//Intrepid_HGRAD_TET_C2_FEM.hpp
81//Intrepid_HGRAD_TET_Cn_FEM.hpp
82//Intrepid_HGRAD_TET_Cn_FEM_ORTH.hpp
83//Intrepid_HGRAD_TET_COMP12_FEM.hpp
84//Intrepid_HGRAD_TRI_C1_FEM.hpp
85//Intrepid_HGRAD_TRI_C2_FEM.hpp
86//Intrepid_HGRAD_TRI_Cn_FEM.hpp
87//Intrepid_HGRAD_TRI_Cn_FEM_ORTH.hpp
88//Intrepid_HGRAD_WEDGE_C1_FEM.hpp
89//Intrepid_HGRAD_WEDGE_C2_FEM.hpp
90//Intrepid_HGRAD_WEDGE_I2_FEM.hpp
91
92// Helper Macro to avoid "unrequested" warnings
93#define MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(level,ename,entry) \
94 {if (level.IsRequested(ename,this) || level.GetKeepFlag(ename,this) != 0) this->Set(level,ename,entry);}
95
96
97
98namespace MueLu {
99
100
101/*********************************************************************************************************/
102namespace MueLuIntrepid {
103inline std::string tolower(const std::string & str) {
104 std::string data(str);
105 std::transform(data.begin(), data.end(), data.begin(), ::tolower);
106 return data;
107}
108
109
110/*********************************************************************************************************/
111 template<class Basis, class LOFieldContainer, class LocalOrdinal, class GlobalOrdinal, class Node>
112 void FindGeometricSeedOrdinals(Teuchos::RCP<Basis> basis, const LOFieldContainer &elementToNodeMap,
113 std::vector<std::vector<LocalOrdinal> > &seeds,
114 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &rowMap,
115 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &columnMap)
116 {
117
118 // For each subcell represented by the elements in elementToNodeMap, we want to identify a globally
119 // unique degree of freedom. Because the other "seed" interfaces in MueLu expect a local ordinal, we
120 // store local ordinals in the resulting seeds container.
121
122 // The approach is as follows. For each element, we iterate through the subcells of the domain topology.
123 // We determine which, if any, of these has the lowest global ID owned that is locally owned. We then insert
124 // the local ID corresponding to this in a vector<set<int>> container whose outer index is the spatial dimension
125 // of the subcell. The set lets us conveniently enforce uniqueness of the stored LIDs.
126
127 shards::CellTopology cellTopo = basis->getBaseCellTopology();
128 int spaceDim = cellTopo.getDimension();
129 seeds.clear();
130 seeds.resize(spaceDim + 1);
131 typedef GlobalOrdinal GO;
132 typedef LocalOrdinal LO;
133
134 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
135 GlobalOrdinal go_invalid = Teuchos::OrdinalTraits<GO>::invalid();
136
137
138 std::vector<std::set<LocalOrdinal>> seedSets(spaceDim+1);
139
140 int numCells = elementToNodeMap.extent(0);
141 auto elementToNodeMap_host = Kokkos::create_mirror_view(elementToNodeMap);
142 Kokkos::deep_copy(elementToNodeMap_host,elementToNodeMap);
143 for (int cellOrdinal=0; cellOrdinal<numCells; cellOrdinal++)
144 {
145 for (int d=0; d<=spaceDim; d++)
146 {
147 int subcellCount = cellTopo.getSubcellCount(d);
148 for (int subcord=0; subcord<subcellCount; subcord++)
149 {
150 int dofCount = basis->getDofCount(d,subcord);
151 if (dofCount == 0) continue;
152 // otherwise, we want to insert the LID corresponding to the least globalID that is locally owned
153 GO leastGlobalDofOrdinal = go_invalid;
154 LO LID_leastGlobalDofOrdinal = lo_invalid;
155 for (int basisOrdinalOrdinal=0; basisOrdinalOrdinal<dofCount; basisOrdinalOrdinal++)
156 {
157 int basisOrdinal = basis->getDofOrdinal(d,subcord,basisOrdinalOrdinal);
158 int colLID = elementToNodeMap_host(cellOrdinal,basisOrdinal);
159 if (colLID != Teuchos::OrdinalTraits<LO>::invalid())
160 {
161 GlobalOrdinal colGID = columnMap.getGlobalElement(colLID);
162 LocalOrdinal rowLID = rowMap.getLocalElement(colGID);
163 if (rowLID != lo_invalid)
164 {
165 if ((leastGlobalDofOrdinal == go_invalid) || (colGID < leastGlobalDofOrdinal))
166 {
167 // replace with rowLID
168 leastGlobalDofOrdinal = colGID;
169 LID_leastGlobalDofOrdinal = rowLID;
170 }
171 }
172 }
173 }
174 if (leastGlobalDofOrdinal != go_invalid)
175 {
176 seedSets[d].insert(LID_leastGlobalDofOrdinal);
177 }
178 }
179 }
180 }
181 for (int d=0; d<=spaceDim;d++)
182 {
183 seeds[d] = std::vector<LocalOrdinal>(seedSets[d].begin(),seedSets[d].end());
184 }
185 }
186
187/*********************************************************************************************************/
188// Syntax [HGRAD|HCURL|HDIV][_| ][HEX|LINE|POLY|PYR|QUAD|TET|TRI|WEDGE][_| ][C|I][1|2|n]
189// Inputs:
190// name - name of the intrepid basis to generate
191// Outputs:
192// degree - order of resulting discretization
193// return value - Intrepid2 basis correspionding to the name
194template<class Scalar,class KokkosExecutionSpace>
195Teuchos::RCP< Intrepid2::Basis<KokkosExecutionSpace,Scalar,Scalar> > BasisFactory(const std::string & name, int & degree)
196{
197 using std::string;
198 using Teuchos::rcp;
199 string myerror("IntrepidBasisFactory: cannot parse string name '"+name+"'");
200
201 // Syntax [HGRAD|HCURL|HDIV][_| ][HEX|LINE|POLY|PYR|QUAD|TET|TRI|WEDGE][_| ][C|I][1|2|n]
202
203 // Get the derivative type
204 size_t pos1 = name.find_first_of(" _");
205 if(pos1==0) throw std::runtime_error(myerror);
206 string deriv = tolower(name.substr(0,pos1));
207 if(deriv!="hgrad" && deriv!="hcurl" && deriv!="hdiv") throw std::runtime_error(myerror);
208
209 // Get the element type
210 pos1++;
211 size_t pos2 = name.find_first_of(" _",pos1);
212 if(pos2==0) throw std::runtime_error(myerror);
213 string el = tolower(name.substr(pos1,pos2-pos1));
214 if(el!="hex" && el!="line" && el!="poly" && el!="pyr" && el!="quad" && el!="tet" && el!="tri" && el!="wedge") throw std::runtime_error(myerror);
215
216 // Get the polynomial type
217 pos2++;
218 string poly = tolower(name.substr(pos2,1));
219 if(poly!="c" && poly!="i") throw std::runtime_error(myerror);
220
221 // Get the degree
222 pos2++;
223 degree=std::stoi(name.substr(pos2,1));
224 if(degree<=0) throw std::runtime_error(myerror);
225
226 // FIXME LATER: Allow for alternative point types for Kirby elements
227 if(deriv=="hgrad" && el=="quad" && poly=="c"){
228 if(degree==1) return rcp(new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<KokkosExecutionSpace,Scalar,Scalar>());
229 else return rcp(new Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar>(degree,Intrepid2::POINTTYPE_EQUISPACED));
230 }
231 else if(deriv=="hgrad" && el=="line" && poly=="c"){
232 if(degree==1) return rcp(new Intrepid2::Basis_HGRAD_LINE_C1_FEM<KokkosExecutionSpace,Scalar,Scalar>());
233 else return rcp(new Intrepid2::Basis_HGRAD_LINE_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar>(degree,Intrepid2::POINTTYPE_EQUISPACED));
234 }
235
236 // Error out
237 throw std::runtime_error(myerror);
238 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
239}
240
241/*********************************************************************************************************/
242// Gets the "lo" nodes nested into a "hi" basis. Only works on quads and lines for a lo basis of p=1
243// Inputs:
244// hi_basis - Higher order Basis
245// Outputs:
246// lo_node_in_hi - std::vector<size_t> of size lo dofs in the reference element, which describes the coindcident hi dots
247// hi_DofCoords - FC<Scalar> of size (#hi dofs, dim) with the coordinate locations of the hi dofs on the reference element
248template<class Scalar,class KokkosDeviceType>
249void IntrepidGetP1NodeInHi(const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space,Scalar,Scalar> >&hi_basis,
250 std::vector<size_t> & lo_node_in_hi,
251 Kokkos::DynRankView<Scalar,KokkosDeviceType> & hi_DofCoords) {
252
253 typedef typename KokkosDeviceType::execution_space KokkosExecutionSpace;
254 // Figure out which unknowns in hi_basis correspond to nodes on lo_basis. This varies by element type.
255 size_t degree = hi_basis->getDegree();
256 lo_node_in_hi.resize(0);
257
258 if(!rcp_dynamic_cast<Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar> >(hi_basis).is_null()) {
259 // HGRAD QUAD Cn: Numbering as per the Kirby convention (straight across, bottom to top)
260 lo_node_in_hi.insert(lo_node_in_hi.end(),{0,degree, (degree+1)*(degree+1)-1, degree*(degree+1)});
261 }
262 else if(!rcp_dynamic_cast<Intrepid2::Basis_HGRAD_LINE_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar> >(hi_basis).is_null()) {
263 // HGRAD LINE Cn: Numbering as per the Kirby convention (straight across)
264 lo_node_in_hi.insert(lo_node_in_hi.end(),{0,degree});
265 }
266 else
267 throw std::runtime_error("IntrepidPCoarsenFactory: Unknown element type");
268
269 // Get coordinates of the hi_basis dof's
270 Kokkos::resize(hi_DofCoords,hi_basis->getCardinality(),hi_basis->getBaseCellTopology().getDimension());
271 hi_basis->getDofCoords(hi_DofCoords);
272}
273
274
275/*********************************************************************************************************/
276// Given a list of candidates picks a definitive list of "representative" higher order nodes for each lo order node via the "smallest GID" rule
277// Input:
278// representative_node_candidates - std::vector<std::vector<size_t> > of lists of "representative candidate" hi dofs for each lo dof
279// hi_elemToNode - FC<LO> containing the high order element-to-node map
280// hi_columnMap - Column map of the higher order matrix
281// Output:
282// lo_elemToHiRepresentativeNode - FC<LO> of size (# elements, # lo dofs per element) listing the hi unknown chosen as the single representative for each lo unknown for counting purposes
283template<class LocalOrdinal, class GlobalOrdinal, class Node, class LOFieldContainer>
284void GenerateLoNodeInHiViaGIDs(const std::vector<std::vector<size_t> > & candidates,const LOFieldContainer & hi_elemToNode,
285 RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & hi_columnMap,
286 LOFieldContainer & lo_elemToHiRepresentativeNode) {
287 typedef GlobalOrdinal GO;
288
289 // Given: A set of "candidate" hi-DOFs to serve as the "representative" DOF for each lo-DOF on the reference element.
290 // Algorithm: For each element, we choose the lowest GID of the candidates for each DOF to generate the lo_elemToHiRepresentativeNode map
291
292 size_t numElem = hi_elemToNode.extent(0);
293 size_t lo_nperel = candidates.size();
294 Kokkos::resize(lo_elemToHiRepresentativeNode,numElem, lo_nperel);
295
296 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
297 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
298 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
299 for(size_t i=0; i<numElem; i++)
300 for(size_t j=0; j<lo_nperel; j++) {
301 if(candidates[j].size() == 1)
302 lo_elemToHiRepresentativeNode_host(i,j) = hi_elemToNode_host(i,candidates[j][0]);
303 else {
304 // First we get the GIDs for each candidate
305 std::vector<GO> GID(candidates[j].size());
306 for(size_t k=0; k<(size_t)candidates[j].size(); k++)
307 GID[k] = hi_columnMap->getGlobalElement(hi_elemToNode_host(i,candidates[j][k]));
308
309 // Find the one with smallest GID
310 size_t which = std::distance(GID.begin(),std::min_element(GID.begin(),GID.end()));
311
312 // Record this
313 lo_elemToHiRepresentativeNode_host(i,j) = hi_elemToNode_host(i,candidates[j][which]);
314 }
315 }
316 Kokkos::deep_copy(lo_elemToHiRepresentativeNode, lo_elemToHiRepresentativeNode_host);
317}
318
319/*********************************************************************************************************/
320// Inputs:
321// hi_elemToNode - FC<LO> containing the high order element-to-node map
322// hi_nodeIsOwned - std::vector<bool> of size hi's column map, which described hi node ownership
323// lo_elemToHiRepresentativeNode - FC<LO> of size (# elements, # lo dofs per element) listing the hi unknown chosen as the single representative for each lo unknown for counting purposes
324// Outputs:
325// lo_elemToNode - FC<LO> containing the low order element-to-node map.
326// lo_nodeIsOwned - std::vector<bool> of size lo's (future) column map, which described lo node ownership
327// hi_to_lo_map - std::vector<LO> of size equal to hi's column map, which contains the lo id each hi idea maps to (or invalid if it doesn't)
328// lo_numOwnedNodes- Number of lo owned nodes
329template <class LocalOrdinal, class LOFieldContainer>
330void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer & hi_elemToNode,
331 const std::vector<bool> & hi_nodeIsOwned,
332 const LOFieldContainer & lo_elemToHiRepresentativeNode,
333 LOFieldContainer & lo_elemToNode,
334 std::vector<bool> & lo_nodeIsOwned,
335 std::vector<LocalOrdinal> & hi_to_lo_map,
336 int & lo_numOwnedNodes) {
337 typedef LocalOrdinal LO;
338 using Teuchos::RCP;
339 // printf("CMS:BuildLoElemToNodeViaRepresentatives: hi_elemToNode.rank() = %d hi_elemToNode.size() = %d\n",hi_elemToNode.rank(), hi_elemToNode.size());
340 size_t numElem = hi_elemToNode.extent(0);
341 size_t hi_numNodes = hi_nodeIsOwned.size();
342 size_t lo_nperel = lo_elemToHiRepresentativeNode.extent(1);
343 Kokkos::resize(lo_elemToNode,numElem, lo_nperel);
344
345 // Start by flagginc the representative nodes
346 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
347 Kokkos::deep_copy(lo_elemToHiRepresentativeNode_host, lo_elemToHiRepresentativeNode);
348 std::vector<bool> is_low_order(hi_numNodes,false);
349 for(size_t i=0; i<numElem; i++)
350 for(size_t j=0; j<lo_nperel; j++) {
351 LO id = lo_elemToHiRepresentativeNode_host(i,j);
352 is_low_order[id] = true; // This can overwrite and that is OK.
353 }
354
355 // Count the number of lo owned nodes, generating a local index for lo nodes
356 lo_numOwnedNodes=0;
357 size_t lo_numNodes=0;
358 hi_to_lo_map.resize(hi_numNodes,Teuchos::OrdinalTraits<LO>::invalid());
359
360 for(size_t i=0; i<hi_numNodes; i++)
361 if(is_low_order[i]) {
362 hi_to_lo_map[i] = lo_numNodes;
363 lo_numNodes++;
364 if(hi_nodeIsOwned[i]) lo_numOwnedNodes++;
365 }
366
367 // Flag the owned lo nodes
368 lo_nodeIsOwned.resize(lo_numNodes,false);
369 for(size_t i=0; i<hi_numNodes; i++) {
370 if(is_low_order[i] && hi_nodeIsOwned[i])
371 lo_nodeIsOwned[hi_to_lo_map[i]]=true;
372 }
373
374 // Translate lo_elemToNode to a lo local index
375 auto lo_elemToNode_host = Kokkos::create_mirror_view(lo_elemToNode);
376 for(size_t i=0; i<numElem; i++)
377 for(size_t j=0; j<lo_nperel; j++)
378 lo_elemToNode_host(i,j) = hi_to_lo_map[lo_elemToHiRepresentativeNode_host(i,j)];
379
380
381 // Check for the [E|T]petra column map ordering property, namely LIDs for owned nodes should all appear first.
382 // Since we're injecting from the higher-order mesh, it should be true, but we should add an error check & throw in case.
383 bool map_ordering_test_passed=true;
384 for(size_t i=0; i<lo_numNodes-1; i++)
385 if(!lo_nodeIsOwned[i] && lo_nodeIsOwned[i+1])
386 map_ordering_test_passed=false;
387
388 if(!map_ordering_test_passed)
389 throw std::runtime_error("MueLu::MueLuIntrepid::BuildLoElemToNodeViaRepresentatives failed map ordering test");
390 Kokkos::deep_copy(lo_elemToNode, lo_elemToNode_host);
391
392}
393
394
395/*********************************************************************************************************/
396// Inputs:
397// hi_elemToNode - FC<LO> containing the high order element-to-node map
398// hi_nodeIsOwned - std::vector<bool> of size hi's column map, which described hi node ownership
399// lo_node_in_hi - std::vector<size_t> of size lo dofs in the reference element, which describes the coindcident hi dots
400// hi_isDirichlet - ArrayView<int> of size of hi's column map, which has a 1 if the unknown is Dirichlet and a 0 if it isn't.
401// Outputs:
402// lo_elemToNode - FC<LO> containing the low order element-to-node map.
403// lo_nodeIsOwned - std::vector<bool> of size lo's (future) column map, which described lo node ownership
404// hi_to_lo_map - std::vector<LO> of size equal to hi's column map, which contains the lo id each hi idea maps to (or invalid if it doesn't)
405// lo_numOwnedNodes- Number of lo owned nodes
406template <class LocalOrdinal, class LOFieldContainer>
407void BuildLoElemToNode(const LOFieldContainer & hi_elemToNode,
408 const std::vector<bool> & hi_nodeIsOwned,
409 const std::vector<size_t> & lo_node_in_hi,
410 const Teuchos::ArrayRCP<const int> & hi_isDirichlet,
411 LOFieldContainer & lo_elemToNode,
412 std::vector<bool> & lo_nodeIsOwned,
413 std::vector<LocalOrdinal> & hi_to_lo_map,
414 int & lo_numOwnedNodes) {
415 typedef LocalOrdinal LO;
416 using Teuchos::RCP;
417 LocalOrdinal LOINVALID = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
418 // printf("CMS:BuildLoElemToNode: hi_elemToNode.rank() = %d hi_elemToNode.size() = %d\n",hi_elemToNode.rank(), hi_elemToNode.size());
419
420 size_t numElem = hi_elemToNode.extent(0);
421 size_t hi_numNodes = hi_nodeIsOwned.size();
422
423 size_t lo_nperel = lo_node_in_hi.size();
424 Kokkos::resize(lo_elemToNode,numElem, lo_nperel);
425
426 // Build lo_elemToNode (in the hi local index ordering) and flag owned ones
427 std::vector<bool> is_low_order(hi_numNodes,false);
428 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
429 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
430 auto lo_elemToNode_host = Kokkos::create_mirror_view(lo_elemToNode);
431 for(size_t i=0; i<numElem; i++)
432 for(size_t j=0; j<lo_nperel; j++) {
433 LO lid = hi_elemToNode_host(i,lo_node_in_hi[j]);
434
435 // Remove Dirichlet
436 if(hi_isDirichlet[lid])
437 lo_elemToNode_host(i,j) = LOINVALID;
438 else {
439 lo_elemToNode_host(i,j) = lid;
440 is_low_order[hi_elemToNode_host(i,lo_node_in_hi[j])] = true; // This can overwrite and that is OK.
441 }
442 }
443
444 // Count the number of lo owned nodes, generating a local index for lo nodes
445 lo_numOwnedNodes=0;
446 size_t lo_numNodes=0;
447 hi_to_lo_map.resize(hi_numNodes,Teuchos::OrdinalTraits<LO>::invalid());
448
449 for(size_t i=0; i<hi_numNodes; i++)
450 if(is_low_order[i]) {
451 hi_to_lo_map[i] = lo_numNodes;
452 lo_numNodes++;
453 if(hi_nodeIsOwned[i]) lo_numOwnedNodes++;
454 }
455
456 // Flag the owned lo nodes
457 lo_nodeIsOwned.resize(lo_numNodes,false);
458 for(size_t i=0; i<hi_numNodes; i++) {
459 if(is_low_order[i] && hi_nodeIsOwned[i])
460 lo_nodeIsOwned[hi_to_lo_map[i]]=true;
461 }
462
463 // Translate lo_elemToNode to a lo local index
464 for(size_t i=0; i<numElem; i++)
465 for(size_t j=0; j<lo_nperel; j++) {
466 if(lo_elemToNode_host(i,j) != LOINVALID)
467 lo_elemToNode_host(i,j) = hi_to_lo_map[lo_elemToNode_host(i,j)];
468 }
469 Kokkos::deep_copy(lo_elemToNode, lo_elemToNode_host);
470
471 // Check for the [E|T]petra column map ordering property, namely LIDs for owned nodes should all appear first.
472 // Since we're injecting from the higher-order mesh, it should be true, but we should add an error check & throw in case.
473 bool map_ordering_test_passed=true;
474 for(size_t i=0; i<lo_numNodes-1; i++)
475 if(!lo_nodeIsOwned[i] && lo_nodeIsOwned[i+1])
476 map_ordering_test_passed=false;
477
478 if(!map_ordering_test_passed)
479 throw std::runtime_error("MueLu::MueLuIntrepid::BuildLoElemToNode failed map ordering test");
480
481}
482
483/*********************************************************************************************************/
484// Generates the lo_columnMap
485// Input:
486// hi_importer - Importer from the hi matrix
487// hi_to_lo_map - std::vector<LO> of size equal to hi's column map, which contains the lo id each hi idea maps to (or invalid if it doesn't)
488// lo_DomainMap - Domain map for the lo matrix
489// lo_columnMapLength - Number of local columns in the lo column map
490// Output:
491// lo_columnMap - Column map of the lower order matrix
492 template <class LocalOrdinal, class GlobalOrdinal, class Node>
493 void GenerateColMapFromImport(const Xpetra::Import<LocalOrdinal,GlobalOrdinal, Node> & hi_importer,const std::vector<LocalOrdinal> &hi_to_lo_map,const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> & lo_domainMap, const size_t & lo_columnMapLength, RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & lo_columnMap) {
494 typedef LocalOrdinal LO;
495 typedef GlobalOrdinal GO;
496 typedef Node NO;
497 typedef Xpetra::Map<LO,GO,NO> Map;
498 typedef Xpetra::Vector<GO,LO,GO,NO> GOVector;
499
500 GO go_invalid = Teuchos::OrdinalTraits<GO>::invalid();
501 LO lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
502
503 RCP<const Map> hi_domainMap = hi_importer.getSourceMap();
504 RCP<const Map> hi_columnMap = hi_importer.getTargetMap();
505 // Figure out the GIDs of my non-owned P1 nodes
506 // HOW: We can build a GOVector(domainMap) and fill the values with either invalid() or the P1 domainMap.GID() for that guy.
507 // Then we can use A's importer to get a GOVector(colMap) with that information.
508
509 // NOTE: This assumes rowMap==colMap and [E|T]petra ordering of all the locals first in the colMap
510 RCP<GOVector> dvec = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(hi_domainMap);
511 {
512 ArrayRCP<GO> dvec_data = dvec->getDataNonConst(0);
513 for(size_t i=0; i<hi_domainMap->getLocalNumElements(); i++) {
514 if(hi_to_lo_map[i]!=lo_invalid) dvec_data[i] = lo_domainMap.getGlobalElement(hi_to_lo_map[i]);
515 else dvec_data[i] = go_invalid;
516 }
517 }
518
519
520 RCP<GOVector> cvec = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(hi_columnMap,true);
521 cvec->doImport(*dvec,hi_importer,Xpetra::ADD);
522
523 // Generate the lo_columnMap
524 // HOW: We can use the local hi_to_lo_map from the GID's in cvec to generate the non-contiguous colmap ids.
525 Array<GO> lo_col_data(lo_columnMapLength);
526 {
527 ArrayRCP<GO> cvec_data = cvec->getDataNonConst(0);
528 for(size_t i=0,idx=0; i<hi_columnMap->getLocalNumElements(); i++) {
529 if(hi_to_lo_map[i]!=lo_invalid) {
530 lo_col_data[idx] = cvec_data[i];
531 idx++;
532 }
533 }
534 }
535
536 lo_columnMap = Xpetra::MapFactory<LO,GO,NO>::Build(lo_domainMap.lib(),Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),lo_col_data(),lo_domainMap.getIndexBase(),lo_domainMap.getComm());
537}
538
539/*********************************************************************************************************/
540// Generates a list of "representative candidate" hi dofs for each lo dof on the reference element. This is to be used in global numbering.
541// Input:
542// basis - The low order basis
543// ReferenceNodeLocations - FC<Scalar> of size (#hidofs, dim) Locations of higher order nodes on the reference element
544// threshold - tolerance for equivalance testing
545// Output:
546// representative_node_candidates - std::vector<std::vector<size_t> > of lists of "representative candidate" hi dofs for each lo dof
547template<class Basis, class SCFieldContainer>
548void GenerateRepresentativeBasisNodes(const Basis & basis, const SCFieldContainer & ReferenceNodeLocations, const double threshold,std::vector<std::vector<size_t> > & representative_node_candidates) {
549 typedef SCFieldContainer FC;
550 typedef typename FC::data_type SC;
551
552 // Evaluate the linear basis functions at the Pn nodes
553 size_t numFieldsHi = ReferenceNodeLocations.extent(0);
554 // size_t dim = ReferenceNodeLocations.extent(1);
555 size_t numFieldsLo = basis.getCardinality();
556
557 FC LoValues("LoValues",numFieldsLo,numFieldsHi);
558
559 basis.getValues(LoValues, ReferenceNodeLocations , Intrepid2::OPERATOR_VALUE);
560
561 Kokkos::fence(); // for kernel in getValues
562
563#if 0
564 printf("** LoValues[%d,%d] **\n",(int)numFieldsLo,(int)numFieldsHi);
565 for(size_t i=0; i<numFieldsLo; i++) {
566 for(size_t j=0; j<numFieldsHi; j++)
567 printf("%6.4e ",LoValues(i,j));
568 printf("\n");
569 }
570 printf("**************\n");fflush(stdout);
571#endif
572
573 representative_node_candidates.resize(numFieldsLo);
574 auto LoValues_host = Kokkos::create_mirror_view(LoValues);
575 Kokkos::deep_copy(LoValues_host, LoValues);
576 for(size_t i=0; i<numFieldsLo; i++) {
577 // 1st pass: find the max value
578 typename Teuchos::ScalarTraits<SC>::magnitudeType vmax = Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<SC>::magnitudeType>::zero();
579 for(size_t j=0; j<numFieldsHi; j++)
580 vmax = std::max(vmax,Teuchos::ScalarTraits<SC>::magnitude(LoValues_host(i,j)));
581
582 // 2nd pass: Find all values w/i threshhold of target
583 for(size_t j=0; j<numFieldsHi; j++) {
584 if(Teuchos::ScalarTraits<SC>::magnitude(vmax - LoValues_host(i,j)) < threshold*vmax)
585 representative_node_candidates[i].push_back(j);
586 }
587 }
588
589 // Sanity check
590 for(size_t i=0; i<numFieldsLo; i++)
591 if(!representative_node_candidates[i].size())
592 throw std::runtime_error("ERROR: GenerateRepresentativeBasisNodes: No candidates found!");
593
594
595}
596
597
598
599}//end MueLu::MueLuIntrepid namespace
600
601
602/*********************************************************************************************************/
603/*********************************************************************************************************/
604/*********************************************************************************************************/
605template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
607 const std::vector<bool> & hi_nodeIsOwned,
608 const SCFieldContainer & hi_DofCoords,
609 const std::vector<size_t> &lo_node_in_hi,
610 const Basis &lo_basis,
611 const std::vector<LocalOrdinal> & hi_to_lo_map,
612 const Teuchos::RCP<const Map> & lo_colMap,
613 const Teuchos::RCP<const Map> & lo_domainMap,
614 const Teuchos::RCP<const Map> & hi_map,
615 Teuchos::RCP<Matrix>& P) const{
616 typedef SCFieldContainer FC;
617 // Evaluate the linear basis functions at the Pn nodes
618 size_t numFieldsHi = hi_elemToNode.extent(1);
619 size_t numFieldsLo = lo_basis.getCardinality();
620 LocalOrdinal LOINVALID = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
621 FC LoValues_at_HiDofs("LoValues_at_HiDofs",numFieldsLo,numFieldsHi);
622 lo_basis.getValues(LoValues_at_HiDofs, hi_DofCoords, Intrepid2::OPERATOR_VALUE);
623 auto LoValues_at_HiDofs_host = Kokkos::create_mirror_view(LoValues_at_HiDofs);
624 Kokkos::deep_copy(LoValues_at_HiDofs_host, LoValues_at_HiDofs);
625 Kokkos::fence(); // for kernel in getValues
626
627 typedef typename Teuchos::ScalarTraits<SC>::halfPrecision SClo;
628 typedef typename Teuchos::ScalarTraits<SClo>::magnitudeType MT;
629 MT effective_zero = Teuchos::ScalarTraits<MT>::eps();
630
631 // Allocate P
632 P = rcp(new CrsMatrixWrap(hi_map,lo_colMap,numFieldsHi)); //FIXLATER: Need faster fill
633 RCP<CrsMatrix> Pcrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
634
635 // Slow-ish fill
636 size_t Nelem=hi_elemToNode.extent(0);
637 std::vector<bool> touched(hi_map->getLocalNumElements(),false);
638 Teuchos::Array<GO> col_gid(1);
639 Teuchos::Array<SC> val(1);
640 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
641 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
642 for(size_t i=0; i<Nelem; i++) {
643 for(size_t j=0; j<numFieldsHi; j++) {
644 LO row_lid = hi_elemToNode_host(i,j);
645 GO row_gid = hi_map->getGlobalElement(row_lid);
646 if(hi_nodeIsOwned[row_lid] && !touched[row_lid]) {
647 for(size_t k=0; k<numFieldsLo; k++) {
648 // Get the local id in P1's column map
649 LO col_lid = hi_to_lo_map[hi_elemToNode_host(i,lo_node_in_hi[k])];
650 if(col_lid==LOINVALID) continue;
651
652 col_gid[0] = {lo_colMap->getGlobalElement(col_lid)};
653 val[0] = LoValues_at_HiDofs_host(k,j);
654
655 // Skip near-zeros
656 if(Teuchos::ScalarTraits<SC>::magnitude(val[0]) >= effective_zero)
657 P->insertGlobalValues(row_gid,col_gid(),val());
658 }
659 touched[row_lid]=true;
660 }
661 }
662 }
663 P->fillComplete(lo_domainMap,hi_map);
664}
665
666/*********************************************************************************************************/
667template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
669 const std::vector<bool> & hi_nodeIsOwned,
670 const SCFieldContainer & hi_DofCoords,
671 const LOFieldContainer & lo_elemToHiRepresentativeNode,
672 const Basis &lo_basis,
673 const std::vector<LocalOrdinal> & hi_to_lo_map,
674 const Teuchos::RCP<const Map> & lo_colMap,
675 const Teuchos::RCP<const Map> & lo_domainMap,
676 const Teuchos::RCP<const Map> & hi_map,
677 Teuchos::RCP<Matrix>& P) const{
678 typedef SCFieldContainer FC;
679 // Evaluate the linear basis functions at the Pn nodes
680 size_t numFieldsHi = hi_elemToNode.extent(1);
681 size_t numFieldsLo = lo_basis.getCardinality();
682 FC LoValues_at_HiDofs("LoValues_at_HiDofs",numFieldsLo,numFieldsHi);
683 lo_basis.getValues(LoValues_at_HiDofs, hi_DofCoords, Intrepid2::OPERATOR_VALUE);
684 auto LoValues_at_HiDofs_host = Kokkos::create_mirror_view(LoValues_at_HiDofs);
685 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
686 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
687 Kokkos::deep_copy(LoValues_at_HiDofs_host, LoValues_at_HiDofs);
688 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
689 Kokkos::deep_copy(lo_elemToHiRepresentativeNode_host, lo_elemToHiRepresentativeNode);
690 Kokkos::fence(); // for kernel in getValues
691
692 typedef typename Teuchos::ScalarTraits<SC>::halfPrecision SClo;
693 typedef typename Teuchos::ScalarTraits<SClo>::magnitudeType MT;
694 MT effective_zero = Teuchos::ScalarTraits<MT>::eps();
695
696 // Allocate P
697 P = rcp(new CrsMatrixWrap(hi_map,lo_colMap,numFieldsHi)); //FIXLATER: Need faster fill
698 RCP<CrsMatrix> Pcrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
699
700 // Slow-ish fill
701 size_t Nelem=hi_elemToNode.extent(0);
702 std::vector<bool> touched(hi_map->getLocalNumElements(),false);
703 Teuchos::Array<GO> col_gid(1);
704 Teuchos::Array<SC> val(1);
705 for(size_t i=0; i<Nelem; i++) {
706 for(size_t j=0; j<numFieldsHi; j++) {
707 LO row_lid = hi_elemToNode_host(i,j);
708 GO row_gid = hi_map->getGlobalElement(row_lid);
709 if(hi_nodeIsOwned[row_lid] && !touched[row_lid]) {
710 for(size_t k=0; k<numFieldsLo; k++) {
711 // Get the local id in P1's column map
712 LO col_lid = hi_to_lo_map[lo_elemToHiRepresentativeNode_host(i,k)];
713 col_gid[0] = {lo_colMap->getGlobalElement(col_lid)};
714 val[0] = LoValues_at_HiDofs_host(k,j);
715
716 // Skip near-zeros
717 if(Teuchos::ScalarTraits<SC>::magnitude(val[0]) >= effective_zero)
718 P->insertGlobalValues(row_gid,col_gid(),val());
719 }
720 touched[row_lid]=true;
721 }
722 }
723 }
724 P->fillComplete(lo_domainMap,hi_map);
725}
726
727/*********************************************************************************************************/
728 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
730 RCP<ParameterList> validParamList = rcp(new ParameterList());
731
732#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
733 SET_VALID_ENTRY("pcoarsen: hi basis");
734 SET_VALID_ENTRY("pcoarsen: lo basis");
735#undef SET_VALID_ENTRY
736
737 validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A used during the prolongator smoothing process");
738
739 validParamList->set< RCP<const FactoryBase> >("Nullspace", Teuchos::null, "Generating factory of the nullspace");
740 validParamList->set< RCP<const FactoryBase> >("pcoarsen: element to node map", Teuchos::null, "Generating factory of the element to node map");
741 return validParamList;
742 }
743
744/*********************************************************************************************************/
745 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
747 Input(fineLevel, "A");
748 Input(fineLevel, "pcoarsen: element to node map");
749 Input(fineLevel, "Nullspace");
750 }
751
752/*********************************************************************************************************/
753 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
755 return BuildP(fineLevel, coarseLevel);
756 }
757
758/*********************************************************************************************************/
759 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
761 FactoryMonitor m(*this, "P Coarsening", coarseLevel);
762 std::string levelIDs = toString(coarseLevel.GetLevelID());
763 const std::string prefix = "MueLu::IntrepidPCoarsenFactory(" + levelIDs + "): ";
764
765 // NOTE: This is hardwired to double on purpose. See the note below.
766 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCi;
767 typedef Kokkos::DynRankView<double,typename Node::device_type> FC;
768
769 // Level Get
770 RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel, "A");
771 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> >(fineLevel, "Nullspace");
772 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Acrs = dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*A);
773
774
775 if (restrictionMode_) {
776 SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
777 A = Utilities::Transpose(*A, true); // build transpose of A explicitly
778 }
779
780 // Find the Dirichlet rows in A
781 std::vector<LocalOrdinal> A_dirichletRows;
782 Utilities::FindDirichletRows(A,A_dirichletRows);
783
784 // Build final prolongator
785 RCP<Matrix> finalP;
786
787 // Reuse pattern if available
788 RCP<ParameterList> APparams = rcp(new ParameterList);
789 if (coarseLevel.IsAvailable("AP reuse data", this)) {
790 GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
791
792 APparams = coarseLevel.Get< RCP<ParameterList> >("AP reuse data", this);
793
794 if (APparams->isParameter("graph"))
795 finalP = APparams->get< RCP<Matrix> >("graph");
796 }
797 const ParameterList& pL = GetParameterList();
798
799 /*******************/
800 // FIXME LATER: Allow these to be manually specified instead of Intrepid
801 // Get the Intrepid bases
802 // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
803 // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet.
804 int lo_degree, hi_degree;
805 RCP<Basis> hi_basis = MueLuIntrepid::BasisFactory<double,typename Node::device_type::execution_space>(pL.get<std::string>("pcoarsen: hi basis"),hi_degree);
806 RCP<Basis> lo_basis = MueLuIntrepid::BasisFactory<double,typename Node::device_type::execution_space>(pL.get<std::string>("pcoarsen: lo basis"),lo_degree);
807
808 // Useful Output
809 GetOStream(Statistics1) << "P-Coarsening from basis "<<pL.get<std::string>("pcoarsen: hi basis")<<" to "<<pL.get<std::string>("pcoarsen: lo basis") <<std::endl;
810
811 /*******************/
812 // Get the higher-order element-to-node map
813 const Teuchos::RCP<FCi> Pn_elemToNode = Get<Teuchos::RCP<FCi> >(fineLevel,"pcoarsen: element to node map");
814
815 // Calculate node ownership (the quick and dirty way)
816 // NOTE: This exploits two things:
817 // 1) domainMap == rowMap
818 // 2) Standard [e|t]petra ordering (namely the local unknowns are always numbered first).
819 // This routine does not work in general.
820 RCP<const Map> rowMap = A->getRowMap();
821 RCP<const Map> colMap = Acrs.getColMap();
822 RCP<const Map> domainMap = A->getDomainMap();
823 int NumProc = rowMap->getComm()->getSize();
824 assert(rowMap->isSameAs(*domainMap));
825 std::vector<bool> Pn_nodeIsOwned(colMap->getLocalNumElements(),false);
826 LO num_owned_rows = 0;
827 for(size_t i=0; i<rowMap->getLocalNumElements(); i++) {
828 if(rowMap->getGlobalElement(i) == colMap->getGlobalElement(i)) {
829 Pn_nodeIsOwned[i] = true;
830 num_owned_rows++;
831 }
832 }
833
834 // Used in all cases
835 FC hi_DofCoords;
836 Teuchos::RCP<FCi> P1_elemToNode = rcp(new FCi());
837
838 std::vector<bool> P1_nodeIsOwned;
839 int P1_numOwnedNodes;
840 std::vector<LO> hi_to_lo_map;
841
842 // Degree-1 variables
843 std::vector<size_t> lo_node_in_hi;
844
845 // Degree-n variables
846 FCi lo_elemToHiRepresentativeNode;
847
848 // Get Dirichlet unknown information
849 RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> > hi_isDirichletRow, hi_isDirichletCol;
850 Utilities::FindDirichletRowsAndPropagateToCols(A,hi_isDirichletRow, hi_isDirichletCol);
851
852#if 0
853 printf("[%d] isDirichletRow = ",A->getRowMap()->getComm()->getRank());
854 for(size_t i=0;i<hi_isDirichletRow->getMap()->getLocalNumElements(); i++)
855 printf("%d ",hi_isDirichletRow->getData(0)[i]);
856 printf("\n");
857 printf("[%d] isDirichletCol = ",A->getRowMap()->getComm()->getRank());
858 for(size_t i=0;i<hi_isDirichletCol->getMap()->getLocalNumElements(); i++)
859 printf("%d ",hi_isDirichletCol->getData(0)[i]);
860 printf("\n");
861 fflush(stdout);
862#endif
863
864 /*******************/
865 if(lo_degree == 1) {
866 // Get reference coordinates and the lo-to-hi injection list for the reference element
867 MueLuIntrepid::IntrepidGetP1NodeInHi(hi_basis,lo_node_in_hi,hi_DofCoords);
868
869 // Generate lower-order element-to-node map
870 MueLuIntrepid::BuildLoElemToNode(*Pn_elemToNode,Pn_nodeIsOwned,lo_node_in_hi,hi_isDirichletCol->getData(0),*P1_elemToNode,P1_nodeIsOwned,hi_to_lo_map,P1_numOwnedNodes);
871 assert(hi_to_lo_map.size() == colMap->getLocalNumElements());
872 }
873 else {
874 // Get lo-order candidates
875 double threshold = 1e-10;
876 std::vector<std::vector<size_t> > candidates;
877 Kokkos::resize(hi_DofCoords,hi_basis->getCardinality(),hi_basis->getBaseCellTopology().getDimension());
878 hi_basis->getDofCoords(hi_DofCoords);
879
880 MueLu::MueLuIntrepid::GenerateRepresentativeBasisNodes<Basis,FC>(*lo_basis,hi_DofCoords,threshold,candidates);
881
882 // Generate the representative nodes
883 MueLu::MueLuIntrepid::GenerateLoNodeInHiViaGIDs(candidates,*Pn_elemToNode,colMap,lo_elemToHiRepresentativeNode);
884 MueLu::MueLuIntrepid::BuildLoElemToNodeViaRepresentatives(*Pn_elemToNode,Pn_nodeIsOwned,lo_elemToHiRepresentativeNode,*P1_elemToNode,P1_nodeIsOwned,hi_to_lo_map,P1_numOwnedNodes);
885 }
886 MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(coarseLevel,"pcoarsen: element to node map",P1_elemToNode);
887
888 /*******************/
889 // Generate the P1_domainMap
890 // HOW: Since we know how many each proc has, we can use the non-uniform contiguous map constructor to do the work for us
891 RCP<const Map> P1_domainMap = MapFactory::Build(rowMap->lib(),Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),P1_numOwnedNodes,rowMap->getIndexBase(),rowMap->getComm());
892 MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(coarseLevel,"CoarseMap",P1_domainMap);
893
894 // Generate the P1_columnMap
895 RCP<const Map> P1_colMap;
896 if(NumProc==1) P1_colMap = P1_domainMap;
897 else MueLuIntrepid::GenerateColMapFromImport<LO,GO,NO>(*Acrs.getCrsGraph()->getImporter(),hi_to_lo_map,*P1_domainMap,P1_nodeIsOwned.size(),P1_colMap);
898
899 /*******************/
900 // Generate the coarsening
901 if(lo_degree == 1)
902 GenerateLinearCoarsening_pn_kirby_to_p1(*Pn_elemToNode,Pn_nodeIsOwned,hi_DofCoords,lo_node_in_hi,*lo_basis,hi_to_lo_map,P1_colMap,P1_domainMap,A->getRowMap(),finalP);
903 else
904 GenerateLinearCoarsening_pn_kirby_to_pm(*Pn_elemToNode,Pn_nodeIsOwned,hi_DofCoords,lo_elemToHiRepresentativeNode,*lo_basis,hi_to_lo_map,P1_colMap,P1_domainMap,A->getRowMap(),finalP);
905
906 /*******************/
907 // Zero out the Dirichlet rows in P
908 Utilities::ZeroDirichletRows(finalP,A_dirichletRows);
909
910 /*******************/
911 // Build the nullspace
912 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P1_domainMap, fineNullspace->getNumVectors());
913 finalP->apply(*fineNullspace,*coarseNullspace,Teuchos::TRANS);
914 Set(coarseLevel, "Nullspace", coarseNullspace);
915
916 // Level Set
917 if (!restrictionMode_) {
918 // The factory is in prolongation mode
919 Set(coarseLevel, "P", finalP);
920
921 APparams->set("graph", finalP);
922 MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(coarseLevel,"AP reuse data",APparams);
923
924 if (IsPrint(Statistics1)) {
925 RCP<ParameterList> params = rcp(new ParameterList());
926 params->set("printLoadBalancingInfo", true);
927 params->set("printCommInfo", true);
928 GetOStream(Statistics1) << PerfUtils::PrintMatrixInfo(*finalP, "P", params);
929 }
930 } else {
931 // The factory is in restriction mode
932 RCP<Matrix> R;
933 {
934 SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
935 R = Utilities::Transpose(*finalP, true);
936 }
937
938 Set(coarseLevel, "R", R);
939
940 if (IsPrint(Statistics2)) {
941 RCP<ParameterList> params = rcp(new ParameterList());
942 params->set("printLoadBalancingInfo", true);
943 params->set("printCommInfo", true);
945 }
946 }
947
948 } //Build()
949
950} //namespace MueLu
951
952#endif // MUELU_IPCFACTORY_DEF_HPP
953
#define SET_VALID_ENTRY(name)
#define MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(level, ename, entry)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Kokkos::DynRankView< LocalOrdinal, typename Node::device_type > LOFieldContainer
void GenerateLinearCoarsening_pn_kirby_to_p1(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const std::vector< size_t > &lo_node_in_hi, const Basis &lo_Basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
Kokkos::DynRankView< double, typename Node::device_type > SCFieldContainer
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Intrepid2::Basis< typename Node::device_type::execution_space, double, double > Basis
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void GenerateLinearCoarsening_pn_kirby_to_pm(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const LOFieldContainer &lo_elemToHiRepresentativeNode, const Basis &lo_basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, Teuchos::RCP< Xpetra::Vector< int, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &isDirichletRow, Teuchos::RCP< Xpetra::Vector< int, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &isDirichletCol)
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, const std::vector< DefaultLocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, std::vector< DefaultLocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
void IntrepidGetP1NodeInHi(const Teuchos::RCP< Intrepid2::Basis< typename KokkosDeviceType::execution_space, Scalar, Scalar > > &hi_basis, std::vector< size_t > &lo_node_in_hi, Kokkos::DynRankView< Scalar, KokkosDeviceType > &hi_DofCoords)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector< std::vector< size_t > > &representative_node_candidates)
void GenerateColMapFromImport(const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &hi_importer, const std::vector< LocalOrdinal > &hi_to_lo_map, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &lo_domainMap, const size_t &lo_columnMapLength, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &lo_columnMap)
void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const LOFieldContainer &lo_elemToHiRepresentativeNode, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
std::string tolower(const std::string &str)
void BuildLoElemToNode(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const std::vector< size_t > &lo_node_in_hi, const Teuchos::ArrayRCP< const int > &hi_isDirichlet, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
void GenerateLoNodeInHiViaGIDs(const std::vector< std::vector< size_t > > &candidates, const LOFieldContainer &hi_elemToNode, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &hi_columnMap, LOFieldContainer &lo_elemToHiRepresentativeNode)
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int &degree)
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.
@ Runtime0
One-liner description of what is happening.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.