|
void | setLocalSizeIndex () |
| Function to calculate the local size indices of sub-domains after MPI domain decomposition. More...
|
|
void | initializeArrays () |
| Function to initialize the arrays used in multi-grid. More...
|
|
void | copyStaggrDerivs () |
| Function to copy the staggered grid derivatives from the grid class to local arrays. More...
|
|
void | setCoefficients () |
| Function to set the coefficients used for calculating laplacian and in smoothing. More...
|
|
void | setStagBounds () |
| Function to set the RectDomain variables for all future references throughout the poisson solver. More...
|
|
virtual void | coarsen () |
| Function to coarsen the grid down the levels of the V-Cycle. More...
|
|
virtual void | prolong () |
| Function to perform prolongation on the array being solved. More...
|
|
virtual void | computeResidual () |
| Function to compute the residual at the start of each V-Cycle. More...
|
|
virtual void | smooth (const int smoothCount) |
| Function to perform smoothing operation on the input array. More...
|
|
virtual real | computeError (const int normOrder) |
| Function to compute the error at the end of each V-Cycle. More...
|
|
virtual void | solve () |
|
virtual void | imposeBC () |
| Function to impose the boundary conditions of Poisson solver at different levels of the V-cycle. More...
|
|
virtual void | createMGSubArrays () |
| Function to create the MPI sub-array data types necessary to transfer data across sub-domains. More...
|
|
virtual void | updatePads (blitz::Array< blitz::Array< real, 3 >, 1 > &data) |
| Function to update the pad points of the local sub-domains at different levels of the V-cycle. More...
|
|
void | vCycle () |
| Function to perform one loop of V-cycle. More...
|
|
|
int | vLevel |
|
int | maxCount |
|
bool | zeroBC |
|
const grid & | mesh |
|
const parser & | inputParams |
|
blitz::Range | all |
|
blitz::Array< blitz::Array< real, 3 >, 1 > | pressureData |
|
blitz::Array< blitz::Array< real, 3 >, 1 > | residualData |
|
blitz::Array< blitz::Array< real, 3 >, 1 > | tmpDataArray |
|
blitz::Array< blitz::Array< real, 3 >, 1 > | smoothedPres |
|
blitz::Array< blitz::RectDomain< 3 >, 1 > | stagFull |
|
blitz::Array< blitz::RectDomain< 3 >, 1 > | stagCore |
|
blitz::Array< int, 1 > | xEnd |
|
blitz::Array< int, 1 > | yEnd |
|
blitz::Array< int, 1 > | zEnd |
|
blitz::Array< int, 1 > | mgSizeArray |
|
blitz::Array< int, 1 > | strideValues |
|
blitz::TinyVector< int, 3 > | localSizeIndex |
|
blitz::Array< MPI_Request, 1 > | recvRequest |
|
blitz::Array< MPI_Status, 1 > | recvStatus |
|
blitz::Array< real, 1 > | hx |
|
blitz::Array< real, 1 > | hy |
|
blitz::Array< real, 1 > | hz |
|
blitz::Array< real, 1 > | hx2 |
|
blitz::Array< real, 1 > | hz2 |
|
blitz::Array< real, 1 > | hzhx |
|
blitz::Array< real, 1 > | hxhy |
|
blitz::Array< real, 1 > | hyhz |
|
blitz::Array< real, 1 > | hxhyhz |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | xixx |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | xix2 |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | etyy |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | ety2 |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | ztzz |
|
blitz::Array< blitz::Array< real, 1 >, 1 > | ztz2 |
|
blitz::Array< MPI_Datatype, 1 > | xMGArray |
|
blitz::Array< MPI_Datatype, 1 > | yMGArray |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgSendLft |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgSendRgt |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgRecvLft |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgRecvRgt |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgSendFrn |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgSendBak |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgRecvFrn |
|
blitz::Array< blitz::TinyVector< int, 3 >, 1 > | mgRecvBak |
|
The base class poisson and its derived classes multigrid_d2 and multigrid_d3.
The class implements the geometric multi-grid method for solving the Poisson equation on a non-uniform grid across MPI decomposed domains for parallel computations. The data structure used by the class for computing multi-grid V-cycles across sub-domains is a blitz array with very wide overlap. When calculating the finite differences at the sub-domain boundaries, at the coarsest level of the V-cycle, data points from very deep within the neighbouring sub-domains are necessary. This is the reason for using a wide pad, that spans up to the nearest node of the adjacent sub-domain at the coarsest mesh. This increases the memory footprint, but doesn't increase the computational time as only a single finite-difference calculation is being done using the pads at all levels of the V-cycle.
All the necessary functions to perform the V-cycle - prolongation, solving at coarsest mesh, smoothening, etc. are implemented within the poisson class.
poisson::poisson |
( |
const grid & |
mesh, |
|
|
const parser & |
solParam |
|
) |
| |
Constructor of the base poisson class.
The short base constructor of the poisson class merely assigns the const references to the grid and parser
class instances being used in the solver.
Moreover, it resizes and populates a local array of multi-grid sizes as used in the grid class.
An array of strides to be used at different V-cycle levels is also generated and stored.
Finally, the maximum allowable number of iterations for the Jacobi iterative solver being used at the
coarsest mesh is set as \f$ N_{max} = N_x \times N_y \times N_z \f$, where \f$N_x\f$, \f$N_y\f$ and \f$N_z\f$
are the number of grid points in the collocated grid at the local sub-domains along x, y and z directions
respectively.
- Parameters
-
mesh | is a const reference to the global data contained in the grid class |
solParam | is a const reference to the user-set parameters contained in the parser class |
void poisson::copyStaggrDerivs |
( |
| ) |
|
|
protected |
Function to copy the staggered grid derivatives from the grid class to local arrays.
Though the grid derivatives in the grid class can be read and accessed, they cannot be used directly
along with the arrays defined in the poisson class as the local arrays have wide pads on both sides.
Therefore, correspondingly wide arrays for grid derivatives are used, into which the staggered grid
derivatives from the grid class are written and stored.
This function serves this purpose of copying the grid derivatives.
The core, publicly accessible function of poisson to compute the solution for the Poisson equation.
The function calls the V-cycle as many times as set by the user.
Before doing so, the input data is transferred into the data-structures used by the poisson class to
perform restrictions and prolongations without copying.
Finally, the computed solution is transferred back from the internal data-structures back into the
scalar field supplied by the calling function.
- Parameters
-
inFn | is a pointer to the plain scalar field (cell-centered) into which the computed soltuion must be transferred |
rhs | is a const reference to the plain scalar field (cell-centered) which contains the RHS for the Poisson equation to solve |
void poisson::setLocalSizeIndex |
( |
| ) |
|
|
protected |
Function to calculate the local size indices of sub-domains after MPI domain decomposition.
For the multi-grid solver, the number of grid nodes must be \f$ 2^N + 1 \f$ to perform V-Cycles
The domain decomposition is also done in such a way that each sub-domain will also have \f$ 2^M + 1 \f$ points.
As a result, the number of processors in each direction must be a power of 2.
In this case, the the value \f$ M \f$ of \f$ 2^M + 1 \f$ can be computed as \f$ M = N - log_2(np) \f$
where \f$ np \f$ is the number of processors along the direction under consideration.
Function to perform one loop of V-cycle.
The V-cycle of restrictions, prolongations and smoothings are performed within this function.
First the input data contained in \ref pressureData is smoothed, after which the residual is computed and stored
in the \ref residualData array.
The restrictions, smoothing, and prolongations are performed on these two arrays subsequently.