1 /**********************************************************************\
2 ______ __ __ _______ _______ __ __ __ __ __ ___
3 / || | | | | ____|| ____|| | | | | \ | | | |/ /
4 | ,----'| | | | | |__ | |__ | | | | | \| | | ' /
5 | | | | | | | __| | __| | | | | | . ` | | <
6 | `----.| `--' | | | | | | `----.| | | |\ | | . \
7 \______| \______/ |__| |__| |_______||__| |__| \__| |__|\__\
11 cufflink is a library for linking numerical methods based on Nvidia's
12 Compute Unified Device Architecture (CUDA™) C/C++ programming language
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/>.
38 Daniel P. Combest. All rights reserved.
39 Modifications by Dominik Christ, Wikki Ltd.
42 diagonal preconditioned conjugate gradient
43 solver for symmetric Matrices using a CUSP CUDA™ based solver.
45 \**********************************************************************/
47 #include "cudaTypes.H"
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>
59 cuspEquationSystem* ces,
60 cudaSolverPerformance* solverPerf,
65 # include "fillCOOMatrix.H"
68 // #include "cufflink/setGPUStorage.H"
70 if(solverPerf->debugCusp)
72 std::cout << " Using amg preconditioner\n";
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;
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);
90 cusp::multiply(A, X, y);
92 //define the normalization factor
93 ValueType normFactor = 1.0;
95 # include "buildNormFactor.H"
98 cusp::blas::axpby(B, y, r, ValueType(1), ValueType(-1));
101 cusp::multiply(M, r, z);
104 cusp::blas::copy(z, p);
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;
116 normR > (solverPerf->tol)
117 && count < (solverPerf->maxIter)
118 && normR/normR0 >= (solverPerf->relTol)
119 || count < solverPerf->minIter
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);
135 cusp::multiply(M, r, z);
137 ValueType rz_old = rz;
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;
152 // End the krylov solver
155 solverPerf->fRes = normR;
156 solverPerf->nIterations = count;
161 solverPerf->fRes<=solverPerf->tol
162 || solverPerf->fRes/solverPerf->iRes<=solverPerf->relTol
165 solverPerf->converged = true;
169 solverPerf->converged = false;
172 // Pass the solution vector back