This file is indexed.

/usr/include/trilinos/AnasaziMinres.hpp is in libtrilinos-anasazi-dev 12.10.1-3.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
// @HEADER
// ***********************************************************************
//
//                 Anasazi: Block Eigensolvers Package
//                 Copyright (2004) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// This library is free software; you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as
// published by the Free Software Foundation; either version 2.1 of the
// License, or (at your option) any later version.
//
// This library is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
// Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
// USA
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
// @HEADER

/*! \file AnasaziMinres.hpp
 *  An implementation of pseudo-block MINRES capable of working with different operators for each RHS.
*/

#ifndef ANASAZI_MINRES_HPP
#define ANASAZI_MINRES_HPP

#include "AnasaziConfigDefs.hpp"
#include "Teuchos_TimeMonitor.hpp"

namespace Anasazi {
namespace Experimental {

template <class Scalar, class MV, class OP>
class PseudoBlockMinres
{
  typedef Anasazi::MultiVecTraits<Scalar,MV>         MVT;
  typedef Anasazi::OperatorTraits<Scalar,MV,OP>      OPT;
  const Scalar ONE; 
  const Scalar ZERO;

public:
  // Constructor
  PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec = Teuchos::null);

  // Set tolerance and maximum iterations
  void setTol(const std::vector<Scalar>& tolerances) { tolerances_ = tolerances; };
  void setMaxIter(const int maxIt) { maxIt_ = maxIt; };

  // Set solution and RHS
  void setSol(Teuchos::RCP<MV> X) { X_ = X; };
  void setRHS(Teuchos::RCP<const MV> B) { B_ = B; };

  // Solve the linear system
  void solve();

private:
  Teuchos::RCP<OP> A_;
  Teuchos::RCP<OP> Prec_;
  Teuchos::RCP<MV> X_;
  Teuchos::RCP<const MV> B_;
  std::vector<Scalar> tolerances_;
  int maxIt_;

  void symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r);

#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
  Teuchos::RCP<Teuchos::Time> AddTime_, ApplyOpTime_, ApplyPrecTime_, AssignTime_, DotTime_, LockTime_, NormTime_, ScaleTime_, TotalTime_;
#endif
};



template <class Scalar, class MV, class OP>
PseudoBlockMinres<Scalar,MV,OP>::PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec) : 
  ONE(Teuchos::ScalarTraits<Scalar>::one()), 
  ZERO(Teuchos::ScalarTraits<Scalar>::zero()) 
{
#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
  AddTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Add Multivectors"); 
  ApplyOpTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Apply Operator");
  ApplyPrecTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Apply Preconditioner");
  AssignTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Assignment (no locking)");
  DotTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Compute Dot Product");
  LockTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Lock Converged Vectors");
  NormTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Compute Norm");
  ScaleTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Scale Multivector");
  TotalTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Total Time");
#endif

  A_ = A;
  Prec_ = Prec;
  maxIt_ = 0;
}



template <class Scalar, class MV, class OP>
void PseudoBlockMinres<Scalar,MV,OP>::solve()
{
  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
    Teuchos::TimeMonitor outertimer( *TotalTime_ );
  #endif

  // Get number of vectors
  int ncols = MVT::GetNumberVecs(*B_);
  int newNumConverged;

  // Declare some variables
  std::vector<Scalar> alpha(ncols), beta, beta1(ncols), phibar, oldBeta(ncols,ZERO), epsln(ncols,ZERO), cs(ncols,-ONE), sn(ncols,ZERO), dbar(ncols,ZERO), oldeps(ncols), delta(ncols), gbar(ncols), phi(ncols), gamma(ncols), tmpvec(ncols);
  std::vector<int> indicesToRemove, newlyConverged, unconvergedIndices(ncols);
  Teuchos::RCP<MV> V, Y, R1, R2, W, W1, W2, locX, scaleHelper, swapHelper;

  // Allocate space for multivectors
  V = MVT::Clone(*B_, ncols);
  Y = MVT::Clone(*B_, ncols);
  R1 = MVT::Clone(*B_, ncols);
  R2 = MVT::Clone(*B_, ncols);
  W = MVT::Clone(*B_, ncols);
  W1 = MVT::Clone(*B_, ncols);
  W2 = MVT::Clone(*B_, ncols);
  scaleHelper = MVT::Clone(*B_, ncols);

  // Reserve space for arrays
  indicesToRemove.reserve(ncols);
  newlyConverged.reserve(ncols);

  // Initialize array of unconverged indices
  for(int i=0; i<ncols; i++)
    unconvergedIndices[i] = i;

  // Get a local copy of X
  // We want the vectors to remain contiguous even as things converge
  locX = MVT::CloneCopy(*X_);

  // Initialize residuals
  // R1 = B - AX
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ );
    #endif
    OPT::Apply(*A_,*locX,*R1);
  }
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *AddTime_ );
    #endif
    MVT::MvAddMv(ONE, *B_, -ONE, *R1, *R1);
  }

  // R2 = R1
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *AssignTime_ );
    #endif
    MVT::Assign(*R1,*R2);
  }

  // Initialize the W's to 0.
  MVT::MvInit (*W);
  MVT::MvInit (*W2);

  // Y = M\R1 (preconditioned residual)
  if(Prec_ != Teuchos::null)
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ );
    #endif

    OPT::Apply(*Prec_,*R1,*Y);
  }
  else
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *AssignTime_ );
    #endif
    MVT::Assign(*R1,*Y);
  }

  // beta1 = sqrt(Y'*R1)
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *DotTime_ );
    #endif
    MVT::MvDot(*Y,*R1, beta1);

    for(size_t i=0; i<beta1.size(); i++)
      beta1[i] = sqrt(beta1[i]);
  }

  // beta = beta1
  beta = beta1;

  // phibar = beta1
  phibar = beta1;

  /////////////////////////////////////////////////////////////////////////////////////////////////
  // Begin Lanczos iterations
  for(int iter=1; iter <= maxIt_; iter++)
  {
    // Test convergence
    indicesToRemove.clear();
    for(int i=0; i<ncols; i++)
    {
      Scalar relres = phibar[i]/beta1[i];
//      std::cout << "relative residual[" << unconvergedIndices[i] << "] at iteration " << iter << " = " << relres << std::endl;

      // If the vector converged, mark it for termination
      // Make sure to add it to 
      if(relres < tolerances_[i])
      {
        indicesToRemove.push_back(i);
      }
    }

    // Check whether anything converged
    newNumConverged = indicesToRemove.size();
    if(newNumConverged > 0)
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *LockTime_ );
      #endif

      // If something converged, stick the converged vectors in X
      newlyConverged.resize(newNumConverged);
      for(int i=0; i<newNumConverged; i++)
        newlyConverged[i] = unconvergedIndices[indicesToRemove[i]];

      Teuchos::RCP<MV> helperLocX = MVT::CloneCopy(*locX,indicesToRemove);

      MVT::SetBlock(*helperLocX,newlyConverged,*X_);

      // If everything has converged, we are done
      if(newNumConverged == ncols)
        return;

      // Update unconverged indices
      std::vector<int> helperVec(ncols - newNumConverged);

      std::set_difference(unconvergedIndices.begin(), unconvergedIndices.end(), newlyConverged.begin(), newlyConverged.end(), helperVec.begin());
      unconvergedIndices = helperVec;

      // Get indices of things we want to keep
      std::vector<int> thingsToKeep(ncols - newNumConverged);
      helperVec.resize(ncols);
      for(int i=0; i<ncols; i++)
        helperVec[i] = i;
      ncols = ncols - newNumConverged;

      std::set_difference(helperVec.begin(), helperVec.end(), indicesToRemove.begin(), indicesToRemove.end(), thingsToKeep.begin());

      // Shrink the multivectors
      Teuchos::RCP<MV> helperMV;
      helperMV = MVT::CloneCopy(*V,thingsToKeep);
      V = helperMV;
      helperMV = MVT::CloneCopy(*Y,thingsToKeep);
      Y = helperMV;
      helperMV = MVT::CloneCopy(*R1,thingsToKeep);
      R1 = helperMV;
      helperMV = MVT::CloneCopy(*R2,thingsToKeep);
      R2 = helperMV;
      helperMV = MVT::CloneCopy(*W,thingsToKeep);
      W = helperMV;
      helperMV = MVT::CloneCopy(*W1,thingsToKeep);
      W1 = helperMV;
      helperMV = MVT::CloneCopy(*W2,thingsToKeep);
      W2 = helperMV;
      helperMV = MVT::CloneCopy(*locX,thingsToKeep);
      locX = helperMV;
      scaleHelper = MVT::Clone(*V,ncols);

      // Shrink the arrays
      alpha.resize(ncols);
      oldeps.resize(ncols);
      delta.resize(ncols);
      gbar.resize(ncols);
      phi.resize(ncols);
      gamma.resize(ncols);
      tmpvec.resize(ncols);
      std::vector<Scalar> helperVecS(ncols);
      for(int i=0; i<ncols; i++)
        helperVecS[i] = beta[thingsToKeep[i]];
      beta = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = beta1[thingsToKeep[i]];
      beta1 = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = phibar[thingsToKeep[i]];
      phibar = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = oldBeta[thingsToKeep[i]];
      oldBeta = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = epsln[thingsToKeep[i]];
      epsln = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = cs[thingsToKeep[i]];
      cs = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = sn[thingsToKeep[i]];
      sn = helperVecS;
      for(int i=0; i<ncols; i++)
        helperVecS[i] = dbar[thingsToKeep[i]];
      dbar = helperVecS;

      // Tell operator about the removed indices
      A_->removeIndices(indicesToRemove);
    }

    // Normalize previous vector
    // V = Y / beta
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*Y,*V);
    }
    for(int i=0; i<ncols; i++)
      tmpvec[i] = ONE/beta[i];
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale (*V, tmpvec);
    }

    // Apply operator
    // Y = AV
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ );
      #endif
      OPT::Apply(*A_, *V, *Y);
    }

    if(iter > 1)
    {
      // Y = Y - beta/oldBeta R1
      for(int i=0; i<ncols; i++)
        tmpvec[i] = beta[i]/oldBeta[i];
      {
        #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
          Teuchos::TimeMonitor lcltimer( *AssignTime_ );
        #endif
        MVT::Assign(*R1,*scaleHelper);
      }
      {
        #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
          Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
        #endif
        MVT::MvScale(*scaleHelper,tmpvec);
      }
      {
        #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
          Teuchos::TimeMonitor lcltimer( *AddTime_ );
        #endif
        MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y);
      }
    }

    // alpha = V'*Y
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *DotTime_ );
      #endif
      MVT::MvDot(*V, *Y, alpha);
    }

    // Y = Y - alpha/beta R2
    for(int i=0; i<ncols; i++)
      tmpvec[i] = alpha[i]/beta[i];
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*R2,*scaleHelper);
    }    
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale(*scaleHelper,tmpvec);
    }
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AddTime_ );
      #endif
      MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y);
    }

    // R1 = R2
    // R2 = Y
    swapHelper = R1;
    R1 = R2;
    R2 = Y;
    Y = swapHelper;

    // Y = M\R2
    if(Prec_ != Teuchos::null)
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ );
      #endif

      OPT::Apply(*Prec_,*R2,*Y);
    }
    else
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*R2,*Y);
    }

    // Get new beta
    // beta = sqrt(R2'*Y)
    oldBeta = beta;
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *DotTime_ );
      #endif
      MVT::MvDot(*Y,*R2, beta);

      for(int i=0; i<ncols; i++)
        beta[i] = sqrt(beta[i]);
    }

    // Apply previous rotation
    oldeps = epsln;
    for(int i=0; i<ncols; i++)
    {
      delta[i] = cs[i]*dbar[i]    + sn[i]*alpha[i];
      gbar[i]  = sn[i]*dbar[i]    - cs[i]*alpha[i];
      epsln[i] =                    sn[i]*beta[i];
      dbar[i]  =                  - cs[i]*beta[i];
    }

    // Compute the next plane rotation
    for(int i=0; i<ncols; i++)
    {
      symOrtho(gbar[i], beta[i], &cs[i], &sn[i], &gamma[i]);

      phi[i] = cs[i]*phibar[i];
      phibar[i] = sn[i]*phibar[i];
    }

    // w1 = w2
    // w2 = w
    swapHelper = W1;
    W1 = W2;
    W2 = W;
    W = swapHelper;

    // W = (V - oldeps*W1 - delta*W2) / gamma
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*W1,*scaleHelper);
    }    
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale(*scaleHelper,oldeps);
    }
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AddTime_ );
      #endif
      MVT::MvAddMv(ONE, *V, -ONE, *scaleHelper, *W);
    }
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*W2,*scaleHelper);
    }    
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale(*scaleHelper,delta);
    }
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AddTime_ );
      #endif
      MVT::MvAddMv(ONE, *W, -ONE, *scaleHelper, *W);
    }
    for(int i=0; i<ncols; i++)
      tmpvec[i] = ONE/gamma[i];
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale(*W,tmpvec);
    }

    // Update X
    // X = X + phi*W
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AssignTime_ );
      #endif
      MVT::Assign(*W,*scaleHelper);
    }    
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
      #endif
      MVT::MvScale(*scaleHelper,phi);
    }
    {
      #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
        Teuchos::TimeMonitor lcltimer( *AddTime_ );
      #endif
      MVT::MvAddMv(ONE, *locX, ONE, *scaleHelper, *locX);
    }
  }

  // Stick unconverged vectors in X
  {
    #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
      Teuchos::TimeMonitor lcltimer( *AssignTime_ );
    #endif
    MVT::SetBlock(*locX,unconvergedIndices,*X_);
  }
}

template <class Scalar, class MV, class OP>
void PseudoBlockMinres<Scalar,MV,OP>::symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r)
{
  const Scalar absA = std::abs(a);
  const Scalar absB = std::abs(b);
  if ( absB == ZERO ) {
    *s = ZERO;
    *r = absA;
    if ( absA == ZERO )
      *c = ONE;
    else
      *c = a / absA;
  } else if ( absA == ZERO ) {
    *c = ZERO;
    *s = b / absB;
    *r = absB;
  } else if ( absB >= absA ) { // && a!=0 && b!=0
    Scalar tau = a / b;
    if ( b < ZERO )
      *s = -ONE / sqrt( ONE+tau*tau );
    else
      *s =  ONE / sqrt( ONE+tau*tau );
    *c = *s * tau;
    *r = b / *s;
  } else { // (absA > absB) && a!=0 && b!=0
    Scalar tau = b / a;
    if ( a < ZERO )
      *c = -ONE / sqrt( ONE+tau*tau );
    else
      *c =  ONE / sqrt( ONE+tau*tau );
    *s = *c * tau;
    *r = a / *c;
  }
}

}} // End of namespace

#endif