Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / cudaSolvers / cudaCG / cgAmg.cu
blobd6ab2807c6791b4fa71fec142042ba0f19b68fb3
1 /**********************************************************************\
2   ______  __    __   _______  _______  __       __   __   __   __  ___
3  /      ||  |  |  | |   ____||   ____||  |     |  | |  \ |  | |  |/  /
4 |  ,----'|  |  |  | |  |__   |  |__   |  |     |  | |   \|  | |  '  /
5 |  |     |  |  |  | |   __|  |   __|  |  |     |  | |  . `  | |    <
6 |  `----.|  `--'  | |  |     |  |     |  `----.|  | |  |\   | |  .  \
7  \______| \______/  |__|     |__|     |_______||__| |__| \__| |__|\__\
9 Cuda For FOAM Link
11 cufflink is a library for linking numerical methods based on Nvidia's
12 Compute Unified Device Architecture (CUDA™) C/C++ programming language
13 and OpenFOAM®.
15 Please note that cufflink is not approved or endorsed by ESI-OpenCFD®
16 Limited, the owner of the OpenFOAM® and OpenCFD® trademarks and
17 producer of OpenFOAM® software.
19 The official web-site of OpenCFD® Limited is www.openfoam.com .
21 ------------------------------------------------------------------------
22 This file is part of cufflink.
24     cufflink is free software: you can redistribute it and/or modify
25     it under the terms of the GNU General Public License as published by
26     the Free Software Foundation, either version 3 of the License, or
27     (at your option) any later version.
29     cufflink is distributed in the hope that it will be useful,
30     but WITHOUT ANY WARRANTY; without even the implied warranty of
31     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
32     GNU General Public License for more details.
34     You should have received a copy of the GNU General Public License
35     along with cufflink.  If not, see <http://www.gnu.org/licenses/>.
37     Author
38     Daniel P. Combest.  All rights reserved.
39     Modifications by Dominik Christ, Wikki Ltd.
41     Description
42     diagonal preconditioned conjugate gradient
43     solver for symmetric Matrices using a CUSP CUDA™ based solver.
45 \**********************************************************************/
47 #include "cudaTypes.H"
49 //CUSP Includes
50 #include <cusp/detail/config.h>
51 #include <cusp/verify.h>
52 #include <cusp/precond/aggregation/smoothed_aggregation.h>
53 #include <cusp/coo_matrix.h>
54 #include <cusp/blas.h>
55 #include <cusp/multiply.h>
57 extern "C" void cgAmg
59     cuspEquationSystem* ces,
60     cudaSolverPerformance* solverPerf,
61     ValueType theta
64     // Populate A
65 #   include "fillCOOMatrix.H"
67     // Set storage
68     // #include "cufflink/setGPUStorage.H"
70     if(solverPerf->debugCusp)
71     {
72       std::cout << "   Using amg preconditioner\n";
73     }
75     cusp::precond::aggregation::
76         smoothed_aggregation<IndexType, ValueType, MemorySpace> M(A);
78     // Start the krylov solver
79     assert(A.num_rows == A.num_cols); // sanity check
81     const size_t N = A.num_rows;
83     // allocate workspace
84     cusp::array1d<ValueType,MemorySpace> y(N);
85     cusp::array1d<ValueType,MemorySpace> z(N);
86     cusp::array1d<ValueType,MemorySpace> r(N);
87     cusp::array1d<ValueType,MemorySpace> p(N);
89     // y <- Ax
90     cusp::multiply(A, X, y);
92     //define the normalization factor
93     ValueType normFactor = 1.0;
95 #   include "buildNormFactor.H"
97     // r <- b - A*x
98     cusp::blas::axpby(B, y, r, ValueType(1), ValueType(-1));
100     // z <- M*r
101     cusp::multiply(M, r, z);
103     // p <- z
104     cusp::blas::copy(z, p);
106     // rz = <r^H, z>
107     ValueType rz = cusp::blas::dotc(r, z);
109     ValueType normR = cusp::blas::nrm2(r)/normFactor;
110     ValueType normR0 = normR;//initial residual
111     solverPerf->iRes = normR0;
112     int count = 0;
114     while
115     (
116         normR > (solverPerf->tol)
117      && count < (solverPerf->maxIter)
118      && normR/normR0 >= (solverPerf->relTol)
119      || count < solverPerf->minIter
120     )
121     {
122         // y <- Ap
123         cusp::multiply(A, p, y);
125         // alpha <- <r,z>/<y,p>
126         ValueType alpha =  rz / cusp::blas::dotc(y, p);
128         // x <- x + alpha * p
129         cusp::blas::axpy(p, X, alpha);
131         // r <- r - alpha * y
132         cusp::blas::axpy(y, r, -alpha);
134         // z <- M*r
135         cusp::multiply(M, r, z);
137         ValueType rz_old = rz;
139         // rz = <r^H, z>
140         rz = cusp::blas::dotc(r, z);
142         // beta <- <r_{i+1},r_{i+1}>/<r,r>
143         ValueType beta = rz / rz_old;
145         // p <- r + beta*p should be p <- z + beta*p
146         cusp::blas::axpby(z, p, p, ValueType(1), beta);
148         normR = cusp::blas::nrm2(r)/normFactor;
150         count++;
151     }
152     // End the krylov solver
154     //final residual
155     solverPerf->fRes = normR;
156     solverPerf->nIterations = count;
158     //converged?
159     if
160     (
161         solverPerf->fRes<=solverPerf->tol
162      || solverPerf->fRes/solverPerf->iRes<=solverPerf->relTol
163     )
164     {
165         solverPerf->converged = true;
166     }
167     else
168     {
169         solverPerf->converged = false;
170     }
172     // Pass the solution vector back
173     ces->X = X;