MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
blas1_d.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

real(kind=8) function dasum (n, x, incx)
 
subroutine daxpy (n, da, dx, incx, dy, incy)
 
subroutine dcopy (n, dx, incx, dy, incy)
 
real(kind=8) function ddot (n, dx, incx, dy, incy)
 
real(kind=8) function dnrm2 (n, x, incx)
 
subroutine drot (n, x, incx, y, incy, c, s)
 
subroutine drotg (sa, sb, c, s)
 
subroutine dscal (n, sa, x, incx)
 
subroutine dswap (n, x, incx, y, incy)
 
integer(kind=4) function idamax (n, dx, incx)
 

Function/Subroutine Documentation

◆ dasum()

real ( kind = 8 ) function dasum ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  x,
integer ( kind = 4 )  incx 
)

Definition at line 1 of file blas1_d.f90.

2 
3  !*****************************************************************************80
4  !
5  !! DASUM takes the sum of the absolute values of a vector.
6  !
7  ! Discussion:
8  !
9  ! This routine uses double precision real arithmetic.
10  !
11  ! Licensing:
12  !
13  ! This code is distributed under the GNU LGPL license.
14  !
15  ! Modified:
16  !
17  ! 15 February 2001
18  !
19  ! Author:
20  !
21  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
22  ! David Kincaid, Fred Krogh.
23  ! FORTRAN90 version by John Burkardt.
24  !
25  ! Reference:
26  !
27  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
28  ! LINPACK User's Guide,
29  ! SIAM, 1979,
30  ! ISBN13: 978-0-898711-72-1,
31  ! LC: QA214.L56.
32  !
33  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
34  ! Algorithm 539,
35  ! Basic Linear Algebra Subprograms for Fortran Usage,
36  ! ACM Transactions on Mathematical Software,
37  ! Volume 5, Number 3, September 1979, pages 308-323.
38  !
39  ! Parameters:
40  !
41  ! Input, integer ( kind = 4 ) N, the number of entries in the vector.
42  !
43  ! Input, real ( kind = 8 ) X(*), the vector to be examined.
44  !
45  ! Input, integer ( kind = 4 ) INCX, the increment between successive
46  ! entries of X. INCX must not be negative.
47  !
48  ! Output, real ( kind = 8 ) DASUM, the sum of the absolute values of X.
49  !
50  implicit none
51 
52  real ( kind = 8 ) dasum
53  integer ( kind = 4 ) incx
54  integer ( kind = 4 ) n
55  real ( kind = 8 ) x(*)
56 
57  dasum = sum( abs( x(1:1+(n-1)*incx:incx) ) )
58 
59  return
real(kind=8) function dasum(n, x, incx)
Definition: blas1_d.f90:2

◆ daxpy()

subroutine daxpy ( integer ( kind = 4 )  n,
real ( kind = 8 )  da,
real ( kind = 8 ), dimension(*)  dx,
integer ( kind = 4 )  incx,
real ( kind = 8 ), dimension(*)  dy,
integer ( kind = 4 )  incy 
)

Definition at line 61 of file blas1_d.f90.

62 
63  !*****************************************************************************80
64  !
65  !! DAXPY computes constant times a vector plus a vector.
66  !
67  ! Discussion:
68  !
69  ! This routine uses double precision real arithmetic.
70  !
71  ! This routine uses unrolled loops for increments equal to one.
72  !
73  ! Licensing:
74  !
75  ! This code is distributed under the GNU LGPL license.
76  !
77  ! Modified:
78  !
79  ! 16 May 2005
80  !
81  ! Author:
82  !
83  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
84  ! David Kincaid, Fred Krogh.
85  ! FORTRAN90 version by John Burkardt.
86  !
87  ! Reference:
88  !
89  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
90  ! LINPACK User's Guide,
91  ! SIAM, 1979,
92  ! ISBN13: 978-0-898711-72-1,
93  ! LC: QA214.L56.
94  !
95  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
96  ! Algorithm 539,
97  ! Basic Linear Algebra Subprograms for Fortran Usage,
98  ! ACM Transactions on Mathematical Software,
99  ! Volume 5, Number 3, September 1979, pages 308-323.
100  !
101  ! Parameters:
102  !
103  ! Input, integer ( kind = 4 ) N, the number of elements in DX and DY.
104  !
105  ! Input, real ( kind = 8 ) DA, the multiplier of DX.
106  !
107  ! Input, real ( kind = 8 ) DX(*), the first vector.
108  !
109  ! Input, integer ( kind = 4 ) INCX, the increment between successive
110  ! entries of DX.
111  !
112  ! Input/output, real ( kind = 8 ) DY(*), the second vector.
113  ! On output, DY(*) has been replaced by DY(*) + DA * DX(*).
114  !
115  ! Input, integer ( kind = 4 ) INCY, the increment between successive
116  ! entries of DY.
117  !
118  implicit none
119 
120  real ( kind = 8 ) da
121  real ( kind = 8 ) dx(*)
122  real ( kind = 8 ) dy(*)
123  integer ( kind = 4 ) i
124  integer ( kind = 4 ) incx
125  integer ( kind = 4 ) incy
126  integer ( kind = 4 ) ix
127  integer ( kind = 4 ) iy
128  integer ( kind = 4 ) m
129  integer ( kind = 4 ) n
130 
131  if ( n <= 0 ) then
132  return
133  end if
134 
135  if ( da == 0.0d+00 ) then
136  return
137  end if
138  !
139  ! Code for unequal increments or equal increments
140  ! not equal to 1.
141  !
142  if ( incx /= 1 .or. incy /= 1 ) then
143 
144  if ( 0 <= incx ) then
145  ix = 1
146  else
147  ix = ( - n + 1 ) * incx + 1
148  end if
149 
150  if ( 0 <= incy ) then
151  iy = 1
152  else
153  iy = ( - n + 1 ) * incy + 1
154  end if
155 
156  do i = 1, n
157  dy(iy) = dy(iy) + da * dx(ix)
158  ix = ix + incx
159  iy = iy + incy
160  end do
161  !
162  ! Code for both increments equal to 1.
163  !
164  else
165 
166  m = mod( n, 4 )
167 
168  dy(1:m) = dy(1:m) + da * dx(1:m)
169 
170  do i = m + 1, n, 4
171  dy(i ) = dy(i ) + da * dx(i )
172  dy(i+1) = dy(i+1) + da * dx(i+1)
173  dy(i+2) = dy(i+2) + da * dx(i+2)
174  dy(i+3) = dy(i+3) + da * dx(i+3)
175  end do
176 
177  end if
178 
179  return

◆ dcopy()

subroutine dcopy ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  dx,
integer ( kind = 4 )  incx,
real ( kind = 8 ), dimension(*)  dy,
integer ( kind = 4 )  incy 
)

Definition at line 181 of file blas1_d.f90.

182 
183  !*****************************************************************************80
184  !
185  !! DCOPY copies a vector X to a vector Y.
186  !
187  ! Discussion:
188  !
189  ! This routine uses double precision real arithmetic.
190  !
191  ! The routine uses unrolled loops for increments equal to one.
192  !
193  ! Licensing:
194  !
195  ! This code is distributed under the GNU LGPL license.
196  !
197  ! Modified:
198  !
199  ! 16 May 2005
200  !
201  ! Author:
202  !
203  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
204  ! David Kincaid, Fred Krogh.
205  ! FORTRAN90 version by John Burkardt.
206  !
207  ! Reference:
208  !
209  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
210  ! LINPACK User's Guide,
211  ! SIAM, 1979,
212  ! ISBN13: 978-0-898711-72-1,
213  ! LC: QA214.L56.
214  !
215  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
216  ! Algorithm 539,
217  ! Basic Linear Algebra Subprograms for Fortran Usage,
218  ! ACM Transactions on Mathematical Software,
219  ! Volume 5, Number 3, September 1979, pages 308-323.
220  !
221  ! Parameters:
222  !
223  ! Input, integer ( kind = 4 ) N, the number of elements in DX and DY.
224  !
225  ! Input, real ( kind = 8 ) DX(*), the first vector.
226  !
227  ! Input, integer ( kind = 4 ) INCX, the increment between successive
228  ! entries of DX.
229  !
230  ! Output, real ( kind = 8 ) DY(*), the second vector.
231  !
232  ! Input, integer ( kind = 4 ) INCY, the increment between successive
233  ! entries of DY.
234  !
235  implicit none
236 
237  integer ( kind = 4 ) i
238  integer ( kind = 4 ) incx
239  integer ( kind = 4 ) incy
240  integer ( kind = 4 ) ix
241  integer ( kind = 4 ) iy
242  integer ( kind = 4 ) m
243  integer ( kind = 4 ) n
244  real ( kind = 8 ) dx(*)
245  real ( kind = 8 ) dy(*)
246 
247  if ( n <= 0 ) then
248  return
249  end if
250 
251  if ( incx == 1 .and. incy == 1 ) then
252 
253  m = mod( n, 7 )
254 
255  if ( m /= 0 ) then
256 
257  dy(1:m) = dx(1:m)
258 
259  end if
260 
261  do i = m+1, n, 7
262  dy(i) = dx(i)
263  dy(i + 1) = dx(i + 1)
264  dy(i + 2) = dx(i + 2)
265  dy(i + 3) = dx(i + 3)
266  dy(i + 4) = dx(i + 4)
267  dy(i + 5) = dx(i + 5)
268  dy(i + 6) = dx(i + 6)
269  end do
270 
271  else
272 
273  if ( 0 <= incx ) then
274  ix = 1
275  else
276  ix = ( -n + 1 ) * incx + 1
277  end if
278 
279  if ( 0 <= incy ) then
280  iy = 1
281  else
282  iy = ( -n + 1 ) * incy + 1
283  end if
284 
285  do i = 1, n
286  dy(iy) = dx(ix)
287  ix = ix + incx
288  iy = iy + incy
289  end do
290 
291  end if
292 
293  return

◆ ddot()

real ( kind = 8 ) function ddot ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  dx,
integer ( kind = 4 )  incx,
real ( kind = 8 ), dimension(*)  dy,
integer ( kind = 4 )  incy 
)

Definition at line 295 of file blas1_d.f90.

296 
297  !*****************************************************************************80
298  !
299  !! DDOT forms the dot product of two vectors.
300  !
301  ! Licensing:
302  !
303  ! This code is distributed under the GNU LGPL license.
304  !
305  ! Modified:
306  !
307  ! 06 June 2014
308  !
309  ! Author:
310  !
311  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
312  ! David Kincaid, Fred Krogh.
313  ! FORTRAN90 version by John Burkardt.
314  !
315  ! Reference:
316  !
317  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
318  ! LINPACK User's Guide,
319  ! SIAM, 1979,
320  ! ISBN13: 978-0-898711-72-1,
321  ! LC: QA214.L56.
322  !
323  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
324  ! Algorithm 539,
325  ! Basic Linear Algebra Subprograms for Fortran Usage,
326  ! ACM Transactions on Mathematical Software,
327  ! Volume 5, Number 3, September 1979, pages 308-323.
328  !
329  ! Parameters:
330  !
331  ! Input, integer ( kind = 4 ) N, the number of entries in the vectors.
332  !
333  ! Input, real ( kind = 8 ) DX(*), the first vector.
334  !
335  ! Input, integer ( kind = 4 ) INCX, the increment between successive
336  ! entries in DX.
337  !
338  ! Input, real ( kind = 8 ) DY(*), the second vector.
339  !
340  ! Input, integer ( kind = 4 ) INCY, the increment between successive
341  ! entries in DY.
342  !
343  ! Output, real ( kind = 8 ) DDOT, the sum of the product of the
344  ! corresponding entries of DX and DY.
345  !
346  implicit none
347 
348  real ( kind = 8 ) ddot
349  real ( kind = 8 ) dx(*)
350  real ( kind = 8 ) dy(*)
351  integer ( kind = 4 ) incx
352  integer ( kind = 4 ) incy
353  integer ( kind = 4 ) n
354  integer ( kind = 4 ) x1
355  integer ( kind = 4 ) xi
356  integer ( kind = 4 ) xn
357  integer ( kind = 4 ) y1
358  integer ( kind = 4 ) yi
359  integer ( kind = 4 ) yn
360  !
361  ! Let the FORTRAN90 function DOT_PRODUCT take care of optimization.
362  !
363  if ( 0 < incx ) then
364  x1 = 1
365  xn = 1 + ( n - 1 ) * incx
366  xi = incx
367  else
368  x1 = 1 + ( n - 1 ) * incx
369  xn = 1
370  xi = - incx
371  end if
372 
373  if ( 0 < incy ) then
374  y1 = 1
375  yn = 1 + ( n - 1 ) * incy
376  yi = incy
377  else
378  y1 = 1 + ( n - 1 ) * incy
379  yn = 1
380  yi = - incy
381  end if
382 
383  ddot = dot_product( dx(x1:xn:xi), dy(y1:yn:yi) )
384 
385  return
real(kind=8) function ddot(n, dx, incx, dy, incy)
Definition: blas1_d.f90:296

◆ dnrm2()

real ( kind = 8 ) function dnrm2 ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  x,
integer ( kind = 4 )  incx 
)

Definition at line 387 of file blas1_d.f90.

388 
389  !*****************************************************************************80
390  !
391  !! DNRM2 returns the euclidean norm of a vector.
392  !
393  ! Discussion:
394  !
395  ! This routine uses double precision real arithmetic.
396  !
397  ! DNRM2 ( X ) = sqrt ( X' * X )
398  !
399  ! Licensing:
400  !
401  ! This code is distributed under the GNU LGPL license.
402  !
403  ! Modified:
404  !
405  ! 16 May 2005
406  !
407  ! Author:
408  !
409  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
410  ! David Kincaid, Fred Krogh.
411  ! FORTRAN90 version by John Burkardt.
412  !
413  ! Reference:
414  !
415  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
416  ! LINPACK User's Guide,
417  ! SIAM, 1979,
418  ! ISBN13: 978-0-898711-72-1,
419  ! LC: QA214.L56.
420  !
421  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
422  ! Algorithm 539,
423  ! Basic Linear Algebra Subprograms for Fortran Usage,
424  ! ACM Transactions on Mathematical Software,
425  ! Volume 5, Number 3, September 1979, pages 308-323.
426  !
427  ! Parameters:
428  !
429  ! Input, integer ( kind = 4 ) N, the number of entries in the vector.
430  !
431  ! Input, real ( kind = 8 ) X(*), the vector whose norm is to be computed.
432  !
433  ! Input, integer ( kind = 4 ) INCX, the increment between successive
434  ! entries of X.
435  !
436  ! Output, real ( kind = 8 ) DNRM2, the Euclidean norm of X.
437  !
438  implicit none
439 
440  real ( kind = 8 ) absxi
441  real ( kind = 8 ) dnrm2
442  integer ( kind = 4 ) incx
443  integer ( kind = 4 ) ix
444  integer ( kind = 4 ) n
445  real ( kind = 8 ) norm
446  real ( kind = 8 ) scale
447  real ( kind = 8 ) ssq
448  real ( kind = 8 ) x(*)
449 
450  if ( n < 1 .or. incx < 1 ) then
451 
452  norm = 0.0d+00
453 
454  else if ( n == 1 ) then
455 
456  norm = abs( x(1) )
457 
458  else
459 
460  scale = 0.0d+00
461  ssq = 1.0d+00
462 
463  do ix = 1, 1 + ( n - 1 ) * incx, incx
464  if ( x(ix) /= 0.0d+00 ) then
465  absxi = abs( x(ix) )
466  if ( scale < absxi ) then
467  ssq = 1.0d+00 + ssq * ( scale / absxi )**2
468  scale = absxi
469  else
470  ssq = ssq + ( absxi / scale )**2
471  end if
472  end if
473  end do
474  norm = scale * sqrt( ssq )
475  end if
476 
477  dnrm2 = norm
478 
479  return
real(kind=8) function dnrm2(n, x, incx)
Definition: blas1_d.f90:388

◆ drot()

subroutine drot ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  x,
integer ( kind = 4 )  incx,
real ( kind = 8 ), dimension(*)  y,
integer ( kind = 4 )  incy,
real ( kind = 8 )  c,
real ( kind = 8 )  s 
)

Definition at line 481 of file blas1_d.f90.

482 
483  !*****************************************************************************80
484  !
485  !! DROT applies a plane rotation.
486  !
487  ! Discussion:
488  !
489  ! This routine uses double precision real arithmetic.
490  !
491  ! Licensing:
492  !
493  ! This code is distributed under the GNU LGPL license.
494  !
495  ! Modified:
496  !
497  ! 08 April 1999
498  !
499  ! Author:
500  !
501  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
502  ! David Kincaid, Fred Krogh.
503  ! FORTRAN90 version by John Burkardt.
504  !
505  ! Reference:
506  !
507  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
508  ! LINPACK User's Guide,
509  ! SIAM, 1979,
510  ! ISBN13: 978-0-898711-72-1,
511  ! LC: QA214.L56.
512  !
513  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
514  ! Algorithm 539,
515  ! Basic Linear Algebra Subprograms for Fortran Usage,
516  ! ACM Transactions on Mathematical Software,
517  ! Volume 5, Number 3, September 1979, pages 308-323.
518  !
519  ! Parameters:
520  !
521  ! Input, integer ( kind = 4 ) N, the number of entries in the vectors.
522  !
523  ! Input/output, real ( kind = 8 ) X(*), one of the vectors to be rotated.
524  !
525  ! Input, integer ( kind = 4 ) INCX, the increment between successive
526  ! entries of X.
527  !
528  ! Input/output, real ( kind = 8 ) Y(*), one of the vectors to be rotated.
529  !
530  ! Input, integer ( kind = 4 ) INCY, the increment between successive
531  ! elements of Y.
532  !
533  ! Input, real ( kind = 8 ) C, S, parameters (presumably the cosine and
534  ! sine of some angle) that define a plane rotation.
535  !
536  implicit none
537 
538  real ( kind = 8 ) c
539  integer ( kind = 4 ) i
540  integer ( kind = 4 ) incx
541  integer ( kind = 4 ) incy
542  integer ( kind = 4 ) ix
543  integer ( kind = 4 ) iy
544  integer ( kind = 4 ) n
545  real ( kind = 8 ) s
546  real ( kind = 8 ) stemp
547  real ( kind = 8 ) x(*)
548  real ( kind = 8 ) y(*)
549 
550  if ( n <= 0 ) then
551 
552  else if ( incx == 1 .and. incy == 1 ) then
553 
554  do i = 1, n
555  stemp = c * x(i) + s * y(i)
556  y(i) = c * y(i) - s * x(i)
557  x(i) = stemp
558  end do
559 
560  else
561 
562  if ( 0 <= incx ) then
563  ix = 1
564  else
565  ix = ( - n + 1 ) * incx + 1
566  end if
567 
568  if ( 0 <= incy ) then
569  iy = 1
570  else
571  iy = ( - n + 1 ) * incy + 1
572  end if
573 
574  do i = 1, n
575  stemp = c * x(ix) + s * y(iy)
576  y(iy) = c * y(iy) - s * x(ix)
577  x(ix) = stemp
578  ix = ix + incx
579  iy = iy + incy
580  end do
581 
582  end if
583 
584  return

◆ drotg()

subroutine drotg ( real ( kind = 8 )  sa,
real ( kind = 8 )  sb,
real ( kind = 8 )  c,
real ( kind = 8 )  s 
)

Definition at line 586 of file blas1_d.f90.

587 
588  !*****************************************************************************80
589  !
590  !! DROTG constructs a Givens plane rotation.
591  !
592  ! Discussion:
593  !
594  ! Given values A and B, this routine computes
595  !
596  ! SIGMA = sign ( A ) if abs ( A ) > abs ( B )
597  ! = sign ( B ) if abs ( A ) <= abs ( B );
598  !
599  ! R = SIGMA * ( A * A + B * B );
600  !
601  ! C = A / R if R is not 0
602  ! = 1 if R is 0;
603  !
604  ! S = B / R if R is not 0,
605  ! 0 if R is 0.
606  !
607  ! The computed numbers then satisfy the equation
608  !
609  ! ( C S ) ( A ) = ( R )
610  ! ( -S C ) ( B ) = ( 0 )
611  !
612  ! The routine also computes
613  !
614  ! Z = S if abs ( A ) > abs ( B ),
615  ! = 1 / C if abs ( A ) <= abs ( B ) and C is not 0,
616  ! = 1 if C is 0.
617  !
618  ! The single value Z encodes C and S, and hence the rotation:
619  !
620  ! If Z = 1, set C = 0 and S = 1;
621  ! If abs ( Z ) < 1, set C = sqrt ( 1 - Z * Z ) and S = Z;
622  ! if abs ( Z ) > 1, set C = 1/ Z and S = sqrt ( 1 - C * C );
623  !
624  ! Licensing:
625  !
626  ! This code is distributed under the GNU LGPL license.
627  !
628  ! Modified:
629  !
630  ! 15 May 2006
631  !
632  ! Author:
633  !
634  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
635  ! David Kincaid, Fred Krogh.
636  ! FORTRAN90 version by John Burkardt.
637  !
638  ! Reference:
639  !
640  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
641  ! LINPACK User's Guide,
642  ! SIAM, 1979,
643  ! ISBN13: 978-0-898711-72-1,
644  ! LC: QA214.L56.
645  !
646  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
647  ! Algorithm 539,
648  ! Basic Linear Algebra Subprograms for Fortran Usage,
649  ! ACM Transactions on Mathematical Software,
650  ! Volume 5, Number 3, September 1979, pages 308-323.
651  !
652  ! Parameters:
653  !
654  ! Input/output, real ( kind = 8 ) SA, SB. On input, SA and SB are the values
655  ! A and B. On output, SA is overwritten with R, and SB is
656  ! overwritten with Z.
657  !
658  ! Output, real ( kind = 8 ) C, S, the cosine and sine of the
659  ! Givens rotation.
660  !
661  implicit none
662 
663  real ( kind = 8 ) c
664  real ( kind = 8 ) r
665  real ( kind = 8 ) roe
666  real ( kind = 8 ) s
667  real ( kind = 8 ) sa
668  real ( kind = 8 ) sb
669  real ( kind = 8 ) scale
670  real ( kind = 8 ) z
671 
672  if ( abs( sb ) < abs( sa ) ) then
673  roe = sa
674  else
675  roe = sb
676  end if
677 
678  scale = abs( sa ) + abs( sb )
679 
680  if ( scale == 0.0d+00 ) then
681  c = 1.0d+00
682  s = 0.0d+00
683  r = 0.0d+00
684  else
685  r = scale * sqrt( ( sa / scale )**2 + ( sb / scale )**2 )
686  r = sign( 1.0d+00, roe ) * r
687  c = sa / r
688  s = sb / r
689  end if
690 
691  if ( 0.0d+00 < abs( c ) .and. abs( c ) <= s ) then
692  z = 1.0d+00 / c
693  else
694  z = s
695  end if
696 
697  sa = r
698  sb = z
699 
700  return

◆ dscal()

subroutine dscal ( integer ( kind = 4 )  n,
real ( kind = 8 )  sa,
real ( kind = 8 ), dimension(*)  x,
integer ( kind = 4 )  incx 
)

Definition at line 702 of file blas1_d.f90.

703 
704  !*****************************************************************************80
705  !
706  !! DSCAL scales a vector by a constant.
707  !
708  ! Discussion:
709  !
710  ! This routine uses double precision real arithmetic.
711  !
712  ! Licensing:
713  !
714  ! This code is distributed under the GNU LGPL license.
715  !
716  ! Modified:
717  !
718  ! 08 April 1999
719  !
720  ! Author:
721  !
722  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
723  ! David Kincaid, Fred Krogh.
724  ! FORTRAN90 version by John Burkardt.
725  !
726  ! Reference:
727  !
728  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
729  ! LINPACK User's Guide,
730  ! SIAM, 1979,
731  ! ISBN13: 978-0-898711-72-1,
732  ! LC: QA214.L56.
733  !
734  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
735  ! Algorithm 539,
736  ! Basic Linear Algebra Subprograms for Fortran Usage,
737  ! ACM Transactions on Mathematical Software,
738  ! Volume 5, Number 3, September 1979, pages 308-323.
739  !
740  ! Parameters:
741  !
742  ! Input, integer ( kind = 4 ) N, the number of entries in the vector.
743  !
744  ! Input, real ( kind = 8 ) SA, the multiplier.
745  !
746  ! Input/output, real ( kind = 8 ) X(*), the vector to be scaled.
747  !
748  ! Input, integer ( kind = 4 ) INCX, the increment between successive
749  ! entries of X.
750  !
751  implicit none
752 
753  integer ( kind = 4 ) i
754  integer ( kind = 4 ) incx
755  integer ( kind = 4 ) ix
756  integer ( kind = 4 ) m
757  integer ( kind = 4 ) n
758  real ( kind = 8 ) sa
759  real ( kind = 8 ) x(*)
760 
761  if ( n <= 0 ) then
762 
763  else if ( incx == 1 ) then
764 
765  m = mod( n, 5 )
766 
767  x(1:m) = sa * x(1:m)
768 
769  do i = m+1, n, 5
770  x(i) = sa * x(i)
771  x(i+1) = sa * x(i+1)
772  x(i+2) = sa * x(i+2)
773  x(i+3) = sa * x(i+3)
774  x(i+4) = sa * x(i+4)
775  end do
776 
777  else
778 
779  if ( 0 <= incx ) then
780  ix = 1
781  else
782  ix = ( - n + 1 ) * incx + 1
783  end if
784 
785  do i = 1, n
786  x(ix) = sa * x(ix)
787  ix = ix + incx
788  end do
789 
790  end if
791 
792  return

◆ dswap()

subroutine dswap ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  x,
integer ( kind = 4 )  incx,
real ( kind = 8 ), dimension(*)  y,
integer ( kind = 4 )  incy 
)

Definition at line 794 of file blas1_d.f90.

795 
796  !*****************************************************************************80
797  !
798  !! DSWAP interchanges two vectors.
799  !
800  ! Discussion:
801  !
802  ! This routine uses double precision real arithmetic.
803  !
804  ! Licensing:
805  !
806  ! This code is distributed under the GNU LGPL license.
807  !
808  ! Modified:
809  !
810  ! 08 April 1999
811  !
812  ! Author:
813  !
814  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
815  ! David Kincaid, Fred Krogh.
816  ! FORTRAN90 version by John Burkardt.
817  !
818  ! Reference:
819  !
820  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
821  ! LINPACK User's Guide,
822  ! SIAM, 1979,
823  ! ISBN13: 978-0-898711-72-1,
824  ! LC: QA214.L56.
825  !
826  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
827  ! Algorithm 539,
828  ! Basic Linear Algebra Subprograms for Fortran Usage,
829  ! ACM Transactions on Mathematical Software,
830  ! Volume 5, Number 3, September 1979, pages 308-323.
831  !
832  ! Parameters:
833  !
834  ! Input, integer ( kind = 4 ) N, the number of entries in the vectors.
835  !
836  ! Input/output, real ( kind = 8 ) X(*), one of the vectors to swap.
837  !
838  ! Input, integer ( kind = 4 ) INCX, the increment between successive
839  ! entries of X.
840  !
841  ! Input/output, real ( kind = 8 ) Y(*), one of the vectors to swap.
842  !
843  ! Input, integer ( kind = 4 ) INCY, the increment between successive
844  ! elements of Y.
845  !
846  implicit none
847 
848  integer ( kind = 4 ) i
849  integer ( kind = 4 ) incx
850  integer ( kind = 4 ) incy
851  integer ( kind = 4 ) ix
852  integer ( kind = 4 ) iy
853  integer ( kind = 4 ) m
854  integer ( kind = 4 ) n
855  real ( kind = 8 ) temp
856  real ( kind = 8 ) x(*)
857  real ( kind = 8 ) y(*)
858 
859  if ( n <= 0 ) then
860 
861  else if ( incx == 1 .and. incy == 1 ) then
862 
863  m = mod( n, 3 )
864 
865  do i = 1, m
866  temp = x(i)
867  x(i) = y(i)
868  y(i) = temp
869  end do
870 
871  do i = m + 1, n, 3
872 
873  temp = x(i)
874  x(i) = y(i)
875  y(i) = temp
876 
877  temp = x(i+1)
878  x(i+1) = y(i+1)
879  y(i+1) = temp
880 
881  temp = x(i+2)
882  x(i+2) = y(i+2)
883  y(i+2) = temp
884 
885  end do
886 
887  else
888 
889  if ( 0 <= incx ) then
890  ix = 1
891  else
892  ix = ( - n + 1 ) * incx + 1
893  end if
894 
895  if ( 0 <= incy ) then
896  iy = 1
897  else
898  iy = ( - n + 1 ) * incy + 1
899  end if
900 
901  do i = 1, n
902  temp = x(ix)
903  x(ix) = y(iy)
904  y(iy) = temp
905  ix = ix + incx
906  iy = iy + incy
907  end do
908 
909  end if
910 
911  return

◆ idamax()

integer ( kind = 4 ) function idamax ( integer ( kind = 4 )  n,
real ( kind = 8 ), dimension(*)  dx,
integer ( kind = 4 )  incx 
)

Definition at line 913 of file blas1_d.f90.

914 
915  !*****************************************************************************80
916  !
917  !! IDAMAX indexes the array element of maximum absolute value.
918  !
919  ! Discussion:
920  !
921  ! This routine uses double precision real arithmetic.
922  !
923  ! Licensing:
924  !
925  ! This code is distributed under the GNU LGPL license.
926  !
927  ! Modified:
928  !
929  ! 08 April 1999
930  !
931  ! Author:
932  !
933  ! Original FORTRAN77 version by Charles Lawson, Richard Hanson,
934  ! David Kincaid, Fred Krogh.
935  ! FORTRAN90 version by John Burkardt.
936  !
937  ! Reference:
938  !
939  ! Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart,
940  ! LINPACK User's Guide,
941  ! SIAM, 1979,
942  ! ISBN13: 978-0-898711-72-1,
943  ! LC: QA214.L56.
944  !
945  ! Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh,
946  ! Algorithm 539,
947  ! Basic Linear Algebra Subprograms for Fortran Usage,
948  ! ACM Transactions on Mathematical Software,
949  ! Volume 5, Number 3, September 1979, pages 308-323.
950  !
951  ! Parameters:
952  !
953  ! Input, integer ( kind = 4 ) N, the number of entries in the vector.
954  !
955  ! Input, real ( kind = 8 ) X(*), the vector to be examined.
956  !
957  ! Input, integer ( kind = 4 ) INCX, the increment between successive
958  ! entries of SX.
959  !
960  ! Output, integer ( kind = 4 ) IDAMAX, the index of the element of SX of
961  ! maximum absolute value.
962  !
963  implicit none
964 
965  real ( kind = 8 ) dmax
966  real ( kind = 8 ) dx(*)
967  integer ( kind = 4 ) i
968  integer ( kind = 4 ) idamax
969  integer ( kind = 4 ) incx
970  integer ( kind = 4 ) ix
971  integer ( kind = 4 ) n
972 
973  idamax = 0
974 
975  if ( n < 1 .or. incx <= 0 ) then
976  return
977  end if
978 
979  idamax = 1
980 
981  if ( n == 1 ) then
982  return
983  end if
984 
985  if ( incx == 1 ) then
986 
987  dmax = abs( dx(1) )
988 
989  do i = 2, n
990  if ( dmax < abs( dx(i) ) ) then
991  idamax = i
992  dmax = abs( dx(i) )
993  end if
994  end do
995 
996  else
997 
998  ix = 1
999  dmax = abs( dx(1) )
1000  ix = ix + incx
1001 
1002  do i = 2, n
1003  if ( dmax < abs( dx(ix) ) ) then
1004  idamax = i
1005  dmax = abs( dx(ix) )
1006  end if
1007  ix = ix + incx
1008  end do
1009 
1010  end if
1011 
1012  return
integer(kind=4) function idamax(n, dx, incx)
Definition: blas1_d.f90:914