Remove trailing whitespace systematically
[foam-extend-3.2.git] / src / cudaSolvers / cudaCG / cgDiag.cu
blobdb4ca0e6271476447b179cc00becc6acc8d3925c
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/diagonal.h>
53 #include <cusp/coo_matrix.h>
54 #include <cusp/ell_matrix.h>
55 #include <cusp/blas.h>
56 #include <cusp/multiply.h>
58 extern "C" void cgDiag
60     cuspEquationSystem* ces,
61     cudaSolverPerformance* solverPerf
64     // Populate A
65 #   include "fillCOOMatrix.H"
67     // Set storage
68     // #include "cufflink/setGPUStorage.H"
70     if(solverPerf->debugCusp)
71     {
72         std::cout << "   Using Cusp_Diagonal preconditioner\n";
73     }
75     cusp::precond::diagonal<ValueType, MemorySpace> M(A);
77     // Start the krylov solver
78     assert(A.num_rows == A.num_cols); // sanity check
80     const size_t N = A.num_rows;
82     // allocate workspace
83     cusp::array1d<ValueType,MemorySpace> y(N);
84     cusp::array1d<ValueType,MemorySpace> z(N);
85     cusp::array1d<ValueType,MemorySpace> r(N);
86     cusp::array1d<ValueType,MemorySpace> p(N);
88     // y <- Ax
89     cusp::multiply(A, X, y);
91     //define the normalization factor
92     ValueType normFactor = 1.0;
94     #include "buildNormFactor.H"
96     // r <- b - A*x
97     cusp::blas::axpby(B, y, r, ValueType(1), ValueType(-1));
99     // z <- M*r
100     cusp::multiply(M, r, z);
102     // p <- z
103     cusp::blas::copy(z, p);
105     // rz = <r^H, z>
106     ValueType rz = cusp::blas::dotc(r, z);
108     ValueType normR = cusp::blas::nrm2(r)/normFactor;
109     ValueType normR0 = normR;   // Initial residual
110     solverPerf->iRes = normR0;
111     int count = 0;
113     while
114     (
115         normR > (solverPerf->tol)
116      && count < (solverPerf->maxIter)
117      && normR/normR0 >= (solverPerf->relTol)
118      || count < solverPerf->minIter
119     )
120     {
121         // y <- Ap
122         cusp::multiply(A, p, y);
124         // alpha <- <r,z>/<y,p>
125         ValueType alpha =  rz / cusp::blas::dotc(y, p);
127         // x <- x + alpha * p
128         cusp::blas::axpy(p, X, alpha);
130         // r <- r - alpha * y
131         cusp::blas::axpy(y, r, -alpha);
133         // z <- M*r
134         cusp::multiply(M, r, z);
136         ValueType rz_old = rz;
138         // rz = <r^H, z>
139         rz = cusp::blas::dotc(r, z);
141         // beta <- <r_{i+1},r_{i+1}>/<r,r>
142         ValueType beta = rz / rz_old;
144         // p <- r + beta*p should be p <- z + beta*p
145         cusp::blas::axpby(z, p, p, ValueType(1), beta);
147         normR = cusp::blas::nrm2(r)/normFactor;
149         count++;
150     }
151     // End the krylov solver
153     // Final residual
154     solverPerf->fRes = normR;
155     solverPerf->nIterations = count;
157     // Converged?
158     if
159     (
160         solverPerf->fRes<=solverPerf->tol
161      || solverPerf->fRes/solverPerf->iRes<=solverPerf->relTol
162     )
163     {
164         solverPerf->converged = true;
165     }
166     else
167     {
168         solverPerf->converged=false;
169     }
171     // Pass the solution vector back
172     ces->X = X;