Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GlobalLinSysPETSc.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File GlobalLinSys.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: GlobalLinSys definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 #include "petscis.h"
39 #include "petscversion.h"
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  /**
46  * @class GlobalLinSysPETSc
47  *
48  * Solves a linear system using PETSc.
49  */
50 
51 
53  const GlobalLinSysKey &pKey,
54  const boost::weak_ptr<ExpList> &pExp,
55  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
56  : GlobalLinSys(pKey, pExp, pLocToGloMap)
57  {
58  // Initialise PETSc
59  PetscInitialize(0, NULL, NULL, NULL);
60 
61  // Create matrix
62  MatCreate(PETSC_COMM_WORLD, &m_matrix);
63  }
64 
66  {
67  }
68 
69  /**
70  * @brief Solve linear system using PETSc.
71  *
72  * The general strategy being a PETSc solve is to:
73  *
74  * - Copy values into the PETSc vector #m_b
75  * - Solve the system #m_ksp and place result into #m_x.
76  * - Scatter results back into #m_locVec using #m_ctx scatter object.
77  * - Copy from #m_locVec to output array #pOutput.
78  */
80  const int pNumRows,
81  const Array<OneD,const NekDouble> &pInput,
82  Array<OneD, NekDouble> &pOutput,
83  const AssemblyMapSharedPtr &locToGloMap,
84  const int pNumDir)
85  {
86  const int nHomDofs = pNumRows - pNumDir;
87 
88  // Populate RHS vector from input
89  VecSetValues(m_b, nHomDofs, &m_reorderedMap[0],
90  &pInput[pNumDir], INSERT_VALUES);
91 
92  // Assemble RHS vector
93  VecAssemblyBegin(m_b);
94  VecAssemblyEnd (m_b);
95 
96  // Do system solve
97  KSPSolve(m_ksp, m_b, m_x);
98 
99  // Scatter results to local vector
100  VecScatterBegin(m_ctx, m_x, m_locVec,
101  INSERT_VALUES, SCATTER_FORWARD);
102  VecScatterEnd (m_ctx, m_x, m_locVec,
103  INSERT_VALUES, SCATTER_FORWARD);
104 
105  // Copy results into output vector
106  PetscScalar *tmp;
107  VecGetArray (m_locVec, &tmp);
108  Vmath::Vcopy (nHomDofs, tmp, 1, &pOutput[pNumDir], 1);
109  VecRestoreArray(m_locVec, &tmp);
110  }
111 
112  /**
113  * @brief Set up PETSc local (equivalent to Nektar++ global) and global
114  * (equivalent to universal) scatter maps.
115  *
116  * These maps are used in GlobalLinSysPETSc::v_SolveLinearSystem to
117  * scatter the solution vector back to each process.
118  */
120  {
121  const int nHomDofs = m_reorderedMap.size();
122 
123  // Create local and global numbering systems for vector
124  IS isGlobal, isLocal;
125  ISCreateGeneral(PETSC_COMM_SELF, nHomDofs, &m_reorderedMap[0],
126  PETSC_COPY_VALUES, &isGlobal);
127  ISCreateStride (PETSC_COMM_SELF, nHomDofs, 0, 1, &isLocal);
128 
129  // Create local vector for output
130  VecCreate (PETSC_COMM_SELF, &m_locVec);
131  VecSetSizes (m_locVec, nHomDofs, PETSC_DECIDE);
132  VecSetFromOptions(m_locVec);
133 
134  // Create scatter context
135  VecScatterCreate (m_x, isGlobal, m_locVec, isLocal, &m_ctx);
136 
137  // Clean up
138  ISDestroy(&isGlobal);
139  ISDestroy(&isLocal);
140  }
141 
142  /**
143  * @brief Calculate a reordering of universal IDs for PETSc.
144  *
145  * PETSc requires a unique, contiguous index of all global and universal
146  * degrees of freedom which represents its position inside the
147  * matrix. Presently Gs does not guarantee this, so this routine
148  * constructs a new universal mapping.
149  */
151  const Array<OneD, const int> &glo2uniMap,
152  const Array<OneD, const int> &glo2unique,
153  const AssemblyMapSharedPtr &pLocToGloMap)
154  {
156  = m_expList.lock()->GetSession()->GetComm();
157 
158  const int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
159  const int nHomDofs = glo2uniMap.num_elements() - nDirDofs;
160  const int nProc = vComm->GetSize();
161  const int rank = vComm->GetRank();
162 
163  int n, cnt;
164 
165  // Count number of unique degrees of freedom on each process.
166  m_nLocal = Vmath::Vsum(nHomDofs, glo2unique + nDirDofs, 1);
167  m_reorderedMap.resize(nHomDofs);
168 
169  // Reduce coefficient counts across all processors.
170  Array<OneD, int> localCounts(nProc, 0), localOffset(nProc, 0);
171  localCounts[rank] = nHomDofs;
172  vComm->AllReduce(localCounts, LibUtilities::ReduceSum);
173 
174  for (n = 1; n < nProc; ++n)
175  {
176  localOffset[n] = localOffset[n-1] + localCounts[n-1];
177  }
178 
179  int totHomDofs = Vmath::Vsum(nProc, localCounts, 1);
180  vector<unsigned int> allUniIds(totHomDofs, 0);
181 
182  // Assemble list of universal IDs
183  for (n = 0; n < nHomDofs; ++n)
184  {
185  int gid = n + nDirDofs;
186  allUniIds[n + localOffset[rank]] = glo2uniMap[gid];
187  }
188 
189  // Reduce this across processors so that each process has a list of
190  // all universal IDs.
191  vComm->AllReduce(allUniIds, LibUtilities::ReduceSum);
192  std::sort(allUniIds.begin(), allUniIds.end());
193  map<int,int> uniIdReorder;
194 
195  // Renumber starting from 0.
196  for (cnt = n = 0; n < allUniIds.size(); ++n)
197  {
198  if (uniIdReorder.count(allUniIds[n]) > 0)
199  {
200  continue;
201  }
202 
203  uniIdReorder[allUniIds[n]] = cnt++;
204  }
205 
206  // Populate reordering map.
207  for (n = 0; n < nHomDofs; ++n)
208  {
209  int gid = n + nDirDofs;
210  int uniId = glo2uniMap[gid];
211  ASSERTL0(uniIdReorder.count(uniId) > 0, "Error in ordering");
212  m_reorderedMap[n] = uniIdReorder[uniId];
213  }
214  }
215 
216  /**
217  * @brief Construct PETSc matrix and vector handles.
218  *
219  * @todo Preallocation should be done at this point, since presently
220  * matrix allocation takes a significant amount of time.
221  */
223  {
224  // CREATE VECTORS
225  VecCreate (PETSC_COMM_WORLD, &m_x);
226  VecSetSizes (m_x, m_nLocal, PETSC_DECIDE);
227  VecSetFromOptions(m_x);
228  VecDuplicate (m_x, &m_b);
229 
230  // CREATE MATRICES
231  MatCreate (PETSC_COMM_WORLD, &m_matrix);
232  MatSetType (m_matrix, MATAIJ);
233  MatSetSizes (m_matrix, m_nLocal, m_nLocal,
234  PETSC_DETERMINE, PETSC_DETERMINE);
235  MatSetFromOptions(m_matrix);
236  MatSetUp (m_matrix);
237  }
238 
239  /**
240  * @brief Set up KSP solver object.
241  *
242  * This is reasonably generic setup -- most solver types can be changed
243  * using the .petscrc file.
244  */
246  {
247  KSPCreate(PETSC_COMM_WORLD, &m_ksp);
248  KSPSetTolerances(
249  m_ksp, tolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
250  KSPSetFromOptions(m_ksp);
251 #if PETSC_VERSION_GE(3,5,0)
252  KSPSetOperators(m_ksp, m_matrix, m_matrix);
253 #else
254  KSPSetOperators(m_ksp, m_matrix, m_matrix, SAME_NONZERO_PATTERN);
255 #endif
256  }
257  }
258 }