MorphoGraphX  2.0-1-227
KrylovMethods.hpp
Go to the documentation of this file.
1 //
2 // This file is part of MorphoGraphX - http://www.MorphoGraphX.org
3 // Copyright (C) 2012-2016 Richard S. Smith and collaborators.
4 //
5 // If you use MorphoGraphX in your work, please cite:
6 // http://dx.doi.org/10.7554/eLife.05864
7 //
8 // MorphoGraphX is free software, and is licensed under under the terms of the
9 // GNU General (GPL) Public License version 2.0, http://www.gnu.org/licenses.
10 //
11 #ifndef KRYLOV_METHODS_HPP
12 #define KRYLOV_METHODS_HPP
13 
14 #include <Geometry.hpp>
15 #include <DistObject.hpp>
16 
17 namespace mgx
18 {
19  // Stabilized Biconjugate Gradient method using distributed matrix
20  template<typename DMatrixT, typename DVectorT>
21  int BiCGStab(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
22  {
23  // Create distributed objects for temps, GPU storage only
24  typename DVectorT::DNhbd &nhbd = x.nhbd();
25  typename DVectorT::DVertex xb_o(nhbd), r_o(nhbd), rnxt_o(nhbd), p_o(nhbd), v_o(nhbd);
26  typename DVectorT::DVertex s_o(nhbd), t_o(nhbd), vec_o(nhbd), r0_o(nhbd), tx_o(nhbd);
27 
28  // Create distributed vectors and matrices
29  DVectorT xb(xb_o), r(r_o), rnxt(rnxt_o), p(p_o), v(v_o);
30  DVectorT s(s_o), t(t_o), vec(vec_o), r0(r0_o), tx(tx_o);
31 
32  // calculate residual
33  r0 = A * x;
34  r0 = b - r0;
35 
36  // first search direction
37  p = r0;
38  r = r0;
39 
40  uint iteration = 0, best = 1;
41  double rnorm, minnorm = 0;
42 
43  do {
44  double alpha, alpha_u, alpha_d;
45  double beta, beta_u, beta_d;
46  double omega, omega_u, omega_d;
47 
48  // v_j
49  v = A * p;
50 
51  // alpha_j
52  alpha_u = r * r0;
53  alpha_d = v * r0;
54  if(alpha_d == 0)
55  alpha = 0;
56  else
57  alpha = alpha_u / alpha_d;
58 
59  // s_j
60  vec = v * alpha;
61  s = r - vec;
62 
63  // t_j
64  t = A * s;
65 
66  // omega
67  omega_u = t * s;
68  omega_d = t * t;
69  if(omega_d == 0)
70  omega = 0;
71  else
72  omega = omega_u / omega_d;
73 
74  // x_(j+1)
75  vec = p * alpha;
76  x += vec;
77  vec = s * omega;
78  x += vec;
79 
80  // r_(j+1)
81  vec = t * omega;
82  rnxt = s - vec;
83 
84  // beta_j
85  beta_u = rnxt * r0;
86  beta_d = r * r0;
87  if(omega == 0 or beta_d == 0)
88  beta = 0;
89  else
90  beta = (alpha*beta_u) / (omega*beta_d);
91 
92  // p_(j+1)
93  vec = v * omega;
94  vec = p - vec;
95  vec *= beta;
96  p = rnxt + vec;
97 
98  // change index of residuals
99  r = rnxt;
100 
101  // Grab norm
102  rnorm = norm(r);
103 
104  // Set min norm, 0 if first time through
105  if(iteration == 0 or rnorm <= minnorm) {
106  minnorm = rnorm;
107  tx = x;
108  best = iteration;
109  }
110  iteration++;
111 
112  if(info > 1)
113  std::cout << iteration << " iterations: " << " norm of residual: " << rnorm << std::endl;
114  }
115  while(!isNan(rnorm) and rnorm < maxnorm and rnorm/100.0 < minnorm and rnorm > tolerance and iteration <= maxiter);
116 
117  // Read result vector back from GPU
118  x = tx;
119 
120  if(info > 0)
121  std::cout<<"BICGSTAB iterations:"<< iteration - 1 << " best iteration:" << best << " norm:" << minnorm << std::endl;
122 
123  return iteration - 1;
124  };
125 
126  // Biconjugate Gradient method using distributed matrix
127  template<typename DMatrixT, typename DVectorT>
128  int BiCG(DMatrixT &A, DMatrixT &AT, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
129  {
130  // Create distributed objects for temps, GPU storage only
131  typename DVectorT::DNhbd &nhbd = x.nhbd();
132  typename DVectorT::DVertex r1_o(nhbd), r2_o(nhbd), p1_o(nhbd), p2_o(nhbd);
133  typename DVectorT::DVertex Ap1_o(nhbd), ATp2_o(nhbd), t_o(nhbd);
134 
135  // Create distributed vectors and matrices
136  DVectorT r1(r1_o), r2(r2_o), p1(p1_o), p2(p2_o), Ap1(Ap1_o), ATp2(ATp2_o), t(t_o);
137 
138  // calculate residual
139  x = 0;
140  r1 = b * -1.0;
141  r2 = r1;
142  p1 = b;
143  p2 = p1;
144  double rr = r1 * r2;
145  if(rr < tolerance) {
146  std::cout << "Initial residual below tolerance" << std::endl;
147  return 0;
148  }
149 
150  uint iteration = 1, best = 1;
151  double minnorm = 0;
152  double oldrr;
153 
154  do {
155  Ap1 = A * p1;
156  ATp2 = AT * p2; // should be A transpose
157 
158  double p2Ap1 = p2 * Ap1;
159  if(info > 0 and p2Ap1 == 0)
160  std::cout << "error" << std::endl;
161 
162  double alpha = rr / p2Ap1;
163  oldrr = rr;
164 
165  t = p1 * alpha;
166  x = x + t;
167  t = Ap1 * alpha;
168  r1 = r1 + t;
169  t = ATp2 * alpha;
170  r2 = r2 + t;
171  rr = r1 * r2;
172 
173  double beta = rr/oldrr;
174  p1 = p1 * beta;
175  p1 = p1 - r1;
176  p2 = p2 * beta;
177  p2 = p2 - r2;
178 
179  // Set min norm, 0 if first time through
180  if(minnorm == 0 or rr <= minnorm) {
181  minnorm = rr;
182  best = iteration;
183  }
184 
185  if(info > 1)
186  std::cout << iteration << " iterations: " << " norm of residual: " << rr << std::endl;
187 
188  iteration++;
189  }
190  while(!isNan(rr) and rr < maxnorm and rr/100.0 < minnorm and rr > tolerance and iteration <= maxiter);
191 
192  if(info > 0)
193  std::cout<<"BICG iterations:"<< iteration - 1 << " best iteration:" << best << " norm:" << minnorm << std::endl;
194 
195  return iteration - 1;
196  };
197 
198  // Conjugate Gradient method using distributed matrix
199  template<typename DMatrixT, typename DVectorT>
200  int CG(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
201  {
202  // Create distributed objects for temps, GPU storage only
203  typename DVectorT::DNhbd &nhbd = x.nhbd();
204  typename DVectorT::DVertex r_o(nhbd), p_o(nhbd), ap_o(nhbd), t_o(nhbd);
205 
206  // Create distributed vectors and matrices
207  DVectorT r(r_o), p(p_o), ap(ap_o);
208 
209  // calculate residual
210  r = A * x;
211  r = b - r;
212  p = r;
213  double rr = r * r;
214  double oldrr = rr;
215 
216  uint iteration = 1, best = 1;
217  double minnorm = 0;
218 
219  do {
220  ap = A * p;
221 
222  double alpha = p * ap;
223  alpha = rr / alpha;
224  x = p * alpha + x;
225  r = ap * -alpha + r;
226 
227  oldrr = rr;
228  rr = r * r;
229  double beta = rr/oldrr;
230  p = p * beta + r;
231 
232  // Set min norm, 0 if first time through
233  if(minnorm == 0 or rr <= minnorm) {
234  minnorm = rr;
235  best = iteration;
236  }
237 
238  if(info > 1)
239  std::cout << iteration << " iterations: " << " norm of residual: " << rr << std::endl;
240 
241  iteration++;
242  }
243  while(!isNan(rr) and rr < maxnorm and rr/100.0 < minnorm and rr > tolerance and iteration <= maxiter);
244 
245  if(info > 0)
246  std::cout<<"CG iterations:"<< iteration - 1 << " best iteration:" << best << " norm:" << minnorm << std::endl;
247 
248  return iteration - 1;
249  };
250 
251  // Conjugate Gradient method with Jacobi pre-conditioner using distributed matrix
252  template<typename DMatrixT, typename DVectorT>
253  int CGPreCond(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
254  {
255  // Create distributed objects for temps, GPU storage only
256  typename DVectorT::DNhbd &nhbd = x.nhbd();
257  typename DVectorT::DVertex r_o(nhbd), p_o(nhbd), ap_o(nhbd), t_o(nhbd), z_o(nhbd);
258  typename DMatrixT::DVertex Mv_o(nhbd);
259  typename DMatrixT::DEdge Me_o(nhbd);
260 
261  // Create distributed vectors and matrices
262  DVectorT r(r_o), p(p_o), ap(ap_o), z(z_o);
263  // Create preconditioner
264  DMatrixT M(Mv_o, Me_o);
265  jacobiPreCond(M, A);
266 
267  // calculate residual
268  r = A * x;
269  r = b - r;
270  z = M * r;
271  p = z;
272  double rz = r * z;
273 
274  uint iteration = 1, best = 1;
275  double minnorm = 0;
276 
277  do {
278  ap = A * p;
279  double alpha = p * ap;
280  alpha = rz / alpha;
281  x = p * alpha + x;
282  r = ap * -alpha + r;
283  z = M * r;
284 
285  double oldrz = rz;
286  rz = r * z;
287  double beta = rz/oldrz;
288  p = p * beta + z;
289 
290  // Set min norm, 0 if first time through
291  if(minnorm == 0 or rz <= minnorm) {
292  minnorm = rz;
293  best = iteration;
294  }
295 
296  if(info > 1)
297  std::cout << iteration << " iterations: " << " norm of residual: " << rz << std::endl;
298 
299  iteration++;
300  }
301  while(!isNan(rz) and rz > tolerance and iteration <= maxiter);
302  //while(!isNan(rz) and rz < maxnorm and rz/100.0 < minnorm and rz > tolerance and iteration <= maxiter);
303 
304  if(info > 0)
305  std::cout<<"CG iterations:"<< iteration - 1 << " best iteration:" << best << " norm:" << minnorm << std::endl;
306 
307  return iteration - 1;
308  };
309 
310  // Preconditioned Stabilized Biconjugate Gradient method using distributed matrix
311  // RSS Does not seem to work, where is the bug?
312  template<typename DMatrixT, typename DVectorT>
313  int BiCGStabPreCond(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
314  {
315  // Create distributed objects for temps, GPU storage only
316  typename DVectorT::DNhbd &nhbd = x.nhbd();
317  typename DVectorT::DVertex xb_o(nhbd), r_o(nhbd), rnxt_o(nhbd), p_o(nhbd), v_o(nhbd), pHat_o(nhbd), sHat_o(nhbd);
318  typename DVectorT::DVertex s_o(nhbd), t_o(nhbd), vec_o(nhbd), r0_o(nhbd), tx_o(nhbd);
319 
320  // Create distributed vectors and matrices
321  DVectorT xb(xb_o), r(r_o), rnxt(rnxt_o), p(p_o), v(v_o), pHat(pHat_o), sHat(sHat_o);;
322  DVectorT s(s_o), t(t_o), vec(vec_o), r0(r0_o), tx(tx_o);
323 
324  // Create preconditioner
325  typename DMatrixT::DVertex Mv_o(nhbd);
326  typename DMatrixT::DEdge Me_o(nhbd);
327 
328  DMatrixT M(Mv_o, Me_o);
329  jacobiPreCond(M, A);
330 
331  // calculate residual
332  r0 = A * x;
333  r0 = b - r0;
334 
335  // first search direction
336  p = r0;
337  r = r0;
338 
339  uint iteration = 0, best = 1;
340  double rnorm, minnorm = 0;
341 
342  do {
343  double alpha, alpha_u, alpha_d;
344  double beta, beta_u, beta_d;
345  double omega, omega_u, omega_d;
346 
347  pHat = M * p;
348 
349  // v_j
350  v = A * pHat;
351 
352  // alpha_j
353  alpha_u = r * r0;
354  alpha_d = v * r0;
355  if(alpha_d == 0)
356  alpha = 0;
357  else
358  alpha = alpha_u / alpha_d;
359 
360  // s_j
361  vec = v * alpha;
362  s = r - vec;
363  sHat = M * s;
364 
365  // t_j
366  t = A * sHat;
367 
368  // omega
369  omega_u = t * s;
370  omega_d = t * t;
371  if(omega_d == 0)
372  omega = 0;
373  else
374  omega = omega_u / omega_d;
375 
376  // x_(j+1)
377  vec = pHat * alpha;
378  x += vec;
379  vec = sHat * omega;
380  x += vec;
381 
382  // r_(j+1)
383  vec = t * omega;
384  rnxt = s - vec;
385 
386  // beta_j
387  beta_u = rnxt * r0;
388  beta_d = r * r0;
389  if(omega == 0 or beta_d == 0)
390  beta = 0;
391  else
392  beta = (alpha*beta_u) / (omega*beta_d);
393 
394  // p_(j+1)
395  vec = v * omega;
396  vec = p - vec;
397  vec *= beta;
398  p = rnxt + vec;
399 
400  // change index of residuals
401  r = rnxt;
402 
403  // Grab norm
404  rnorm = norm(r);
405 
406  // Set min norm, 0 if first time through
407  if(iteration == 0 or rnorm <= minnorm) {
408  minnorm = rnorm;
409  tx = x;
410  best = iteration;
411  }
412  iteration++;
413 
414  if(info > 1)
415  std::cout << iteration << " iterations: " << " norm of residual: " << rnorm << std::endl;
416  }
417  while(!isNan(rnorm) and rnorm < maxnorm and rnorm > tolerance and iteration <= maxiter);
418  //while(!isNan(rnorm) and rnorm < maxnorm and rnorm/100.0 < minnorm and rnorm > tolerance and iteration <= maxiter);
419 
420  // Read result vector back from GPU
421  x = tx;
422 
423  if(info > 0)
424  std::cout<<"BICGSTAB iterations:"<< iteration - 1 << " best iteration:" << best << " norm:" << minnorm << std::endl;
425 
426  return iteration - 1;
427  };
428 }
429 #endif
mgx::jacobiPreCond
void jacobiPreCond(DMatrix< DistNhbdT, MatrixT, VectorT, ScalarT > &m, DMatrix< DistNhbdT, MatrixT, VectorT, ScalarT > &a)
Definition: DistMatrix.hpp:557
mgx::uint
unsigned int uint
Definition: Geometry.hpp:41
mgx::BiCGStab
int BiCGStab(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
Definition: KrylovMethods.hpp:21
mgx::BiCG
int BiCG(DMatrixT &A, DMatrixT &AT, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
Definition: KrylovMethods.hpp:128
mgx::CGPreCond
int CGPreCond(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
Definition: KrylovMethods.hpp:253
Geometry.hpp
mgx
Distributed matrix library.
Definition: Assert.hpp:26
mgx::CG
int CG(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
Definition: KrylovMethods.hpp:200
mgx::BiCGStabPreCond
int BiCGStabPreCond(DMatrixT &A, DVectorT &x, DVectorT &b, double tolerance, uint maxiter=500, double maxnorm=1e100, int info=0)
Definition: KrylovMethods.hpp:313
DistObject.hpp
mgx::norm
ScalarT norm(DVector< DistNhbdT, MatrixT, VectorT, ScalarT > &v)
Definition: DistMatrix.hpp:540
mgx::isNan
CU_HOST_DEVICE bool isNan(float s)
Definition: Geometry.hpp:148