Project
Loading...
Searching...
No Matches
PoissonSolver.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
21
22#ifndef ALICEO2_TPC_POISSONSOLVER_H_
23#define ALICEO2_TPC_POISSONSOLVER_H_
24
28#include <vector>
29
30namespace o2
31{
32namespace tpc
33{
34
35template <typename DataT>
36class Vector3D;
37
38template <typename DataT>
39class DataContainer3D;
40
45
47template <typename DataT = double>
49{
50 public:
54
56 PoissonSolver(const RegularGrid& gridProperties) : mGrid3D{gridProperties} {};
57
74 void poissonSolver3D(DataContainer& matricesV, const DataContainer& matricesCharge, const int symmetry);
75
82 void poissonSolver2D(DataContainer& matricesV, const DataContainer& matricesCharge);
83
84 DataT getSpacingZ() const { return mGrid3D.getSpacingZ(); }
85 DataT getSpacingR() const { return mGrid3D.getSpacingR(); }
86 DataT getSpacingPhi() const { return mGrid3D.getSpacingPhi(); }
87
88 static void setConvergenceError(const DataT error) { sConvergenceError = error; }
89
90 static DataT getConvergenceError() { return sConvergenceError; }
91
93 static int getNThreads() { return sNThreads; }
94
96 static void setNThreads(int nThreads) { sNThreads = nThreads; }
97
98 private:
99 const RegularGrid& mGrid3D{};
100 const ParamSpaceCharge mParamGrid{mGrid3D.getParamSC()};
101 inline static DataT sConvergenceError{1e-6};
102 static constexpr DataT INVTWOPI = 1. / o2::constants::math::TwoPI;
103 inline static int sNThreads{4};
104
106 static DataT getGridSizePhiInv();
107
112 DataT getConvergenceError(const Vector& matricesCurrentV, Vector& prevArrayV) const;
113
136 //
139 void poissonMultiGrid3D2D(DataContainer& matricesV, const DataContainer& matricesCharge, const int symmetry);
140
163 void poissonMultiGrid3D(DataContainer& matricesV, const DataContainer& matricesCharge, const int symmetry);
164
173 void poissonMultiGrid2D(DataContainer& matricesV, const DataContainer& matricesCharge, const int iPhi = 0);
174
190 void restrict2D(Vector& matricesCurrentCharge, const Vector& residue, const int tnRRow, const int tnZColumn, const int iphi) const;
191
209 void restrict3D(Vector& matricesCurrentCharge, const Vector& residue, const int tnRRow, const int tnZColumn, const int newPhiSlice, const int oldPhiSlice) const;
210
221 void restrictBoundary3D(Vector& matricesCurrentCharge, const Vector& residue, const int tnRRow, const int tnZColumn, const int newPhiSlice, const int oldPhiSlice) const;
222
241 void relax3D(Vector& matricesCurrentV, const Vector& matricesCurrentCharge, const int tnRRow, const int tnZColumn, const int iPhi, const int symmetry, const DataT h2, const DataT tempRatioZ,
242 const std::vector<DataT>& coefficient1, const std::vector<DataT>& coefficient2, const std::vector<DataT>& coefficient3, const std::vector<DataT>& coefficient4) const;
243
261 void relax2D(Vector& matricesCurrentV, const Vector& matricesCurrentCharge, const int tnRRow, const int tnZColumn, const DataT h2, const DataT tempFourth, const DataT tempRatio,
262 std::vector<DataT>& coefficient1, std::vector<DataT>& coefficient2);
263
279 void interp2D(Vector& matricesCurrentV, const Vector& matricesCurrentVC, const int tnRRow, const int tnZColumn, const int iphi) const;
280
298 void interp3D(Vector& matricesCurrentV, const Vector& matricesCurrentVC, const int tnRRow, const int tnZColumn, const int newPhiSlice, const int oldPhiSlice) const;
299
312 void addInterp3D(Vector& matricesCurrentV, const Vector& matricesCurrentVC, const int tnRRow, const int tnZColumn, const int newPhiSlice, const int oldPhiSlice) const;
313
324 void addInterp2D(Vector& matricesCurrentV, const Vector& matricesCurrentVC, const int tnRRow, const int tnZColumn, const int tnPhi) const;
325
353 void vCycle3D2D(const int symmetry, const int gridFrom, const int gridTo, const int nPre, const int nPost, const DataT ratioZ, const DataT ratioPhi, std::vector<Vector>& tvArrayV,
354 std::vector<Vector>& tvCharge, std::vector<Vector>& tvResidue, std::vector<DataT>& coefficient1, std::vector<DataT>& coefficient2, std::vector<DataT>& coefficient3,
355 std::vector<DataT>& coefficient4, std::vector<DataT>& inverseCoefficient4) const;
356
382 void vCycle3D(const int symmetry, const int gridFrom, const int gridTo, const int nPre, const int nPost, const DataT ratioZ, std::vector<Vector>& tvArrayV, std::vector<Vector>& tvCharge,
383 std::vector<Vector>& tvResidue, std::vector<DataT>& coefficient1, std::vector<DataT>& coefficient2, std::vector<DataT>& coefficient3,
384 std::vector<DataT>& coefficient4, std::vector<DataT>& inverseCoefficient4) const;
385
401 void vCycle2D(const int gridFrom, const int gridTo, const int nPre, const int nPost, const DataT gridSizeR, const DataT ratio, std::vector<Vector>& tvArrayV,
402 std::vector<Vector>& tvCharge, std::vector<Vector>& tvResidue);
403
420 void wCycle2D(const int gridFrom, const int gridTo, const int gamma, const int nPre, const int nPost, const DataT gridSizeR, const DataT ratio,
421 std::vector<Vector>& tvArrayV, std::vector<Vector>& tvCharge, std::vector<Vector>& tvResidue);
422
444 void residue3D(Vector& residue, const Vector& matricesCurrentV, const Vector& matricesCurrentCharge, const int tnRRow, const int tnZColumn, const int tnPhi, const int symmetry, const DataT ih2, const DataT tempRatioZ,
445 const std::vector<DataT>& coefficient1, const std::vector<DataT>& coefficient2, const std::vector<DataT>& coefficient3, const std::vector<DataT>& inverseCoefficient4) const;
446
465 void residue2D(Vector& residue, const Vector& matricesCurrentV, const Vector& matricesCurrentCharge, const int tnRRow, const int tnZColumn, const DataT ih2, const DataT inverseTempFourth,
466 const DataT tempRatio, std::vector<DataT>& coefficient1, std::vector<DataT>& coefficient2);
467
475 void restrictBoundary2D(Vector& matricesCurrentCharge, const Vector& residue, const int tnRRow, const int tnZColumn, const int tnPhi) const;
476
477 // calculate coefficients
478 void calcCoefficients(unsigned int from, unsigned int to, const DataT h, const DataT tempRatioZ, const DataT tempRatioPhi, std::vector<DataT>& coefficient1,
479 std::vector<DataT>& coefficient2, std::vector<DataT>& coefficient3, std::vector<DataT>& coefficient4) const;
480
481 // calculate coefficients for 2D poisson solver
482 void calcCoefficients2D(unsigned int from, unsigned int to, const DataT h, std::vector<DataT>& coefficient1, std::vector<DataT>& coefficient2) const;
483
487 bool isPowerOfTwo(const int i) const
488 {
489 return ((i > 0) && !(i & (i - 1)));
490 };
491};
492
493} // namespace tpc
494} // namespace o2
495
496#endif
int32_t i
useful math constants
Definition of RegularGrid3D class.
Class for time synchronization of RawReader instances.
DataT getSpacingZ() const
DataT getSpacingR() const
static void setConvergenceError(const DataT error)
DataContainer3D< DataT > DataContainer
DataT getSpacingPhi() const
static void setNThreads(int nThreads)
set the number of threads used for some of the calculations
PoissonSolver(const RegularGrid &gridProperties)
default constructor
Vector3D< DataT > Vector
static DataT getConvergenceError()
RegularGrid3D< DataT > RegularGrid
static int getNThreads()
get the number of threads used for some of the calculations
const auto & getParamSC() const
get max index in z direction
DataT getSpacingPhi() const
DataT getSpacingR() const
DataT getSpacingZ() const
this is a simple vector class which is used in the poisson solver class
Definition Vector3D.h:31
constexpr float TwoPI
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< T >, ROOT::Math::DefaultCoordinateSystemTag > Vector3D
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...