/[MITgcm]/MITgcm/lsopt/dgscale.F
ViewVC logotype

Contents of /MITgcm/lsopt/dgscale.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (show annotations) (download)
Fri Nov 15 04:03:24 2002 UTC (21 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint47e_post, checkpoint47c_post, checkpoint50c_post, checkpoint48e_post, checkpoint50c_pre, checkpoint51o_pre, checkpoint51l_post, checkpoint48i_post, checkpoint50d_pre, checkpoint51, checkpoint50, checkpoint52, checkpoint50d_post, checkpoint50b_pre, checkpoint51f_post, checkpoint48b_post, checkpoint51d_post, checkpoint48c_pre, checkpoint47d_pre, checkpoint51t_post, checkpoint51n_post, checkpoint51s_post, checkpoint47a_post, checkpoint48d_pre, checkpoint51j_post, checkpoint47i_post, checkpoint51n_pre, checkpoint47d_post, checkpoint48d_post, checkpoint48f_post, checkpoint51l_pre, checkpoint48h_post, checkpoint51q_post, checkpoint51b_pre, checkpoint47g_post, checkpoint51h_pre, checkpoint48a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint47j_post, branch-exfmods-tag, branchpoint-genmake2, checkpoint51r_post, checkpoint48c_post, checkpoint51i_post, checkpoint51b_post, checkpoint51c_post, checkpoint47b_post, checkpoint50g_post, checkpoint52a_pre, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51i_pre, checkpoint47f_post, checkpoint50e_post, checkpoint51e_post, checkpoint47, checkpoint48, checkpoint49, checkpoint51o_post, checkpoint51f_pre, checkpoint48g_post, checkpoint47h_post, checkpoint52a_post, checkpoint51g_post, ecco_c52_e35, checkpoint50b_post, checkpoint51m_post, checkpoint51a_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-exfmods-curt, branch-genmake2, branch-nonh, tg2-branch, checkpoint51n_branch
Changes since 1.1: +74 -0 lines
o Incorporating QNVS line search routines into MITgcm
  (this is separate code, not compiled with MITgcm,
  and therefore not under pkg)
  - lsopt/
  - optim/

1
2 subroutine dgscale( nn, gold, xdiff, diag, rmin )
3
4
5 c ==================================================================
6 c SUBROUTINE dgscale
7 c ==================================================================
8 c
9 c o computes new preconditioner and writes it to OPWARMD
10 c
11 c o started: ??? not reproducible
12 c
13 c o Version: 2.1.0, 02-Mar-2000: Patrick Heimbach, MIT/EAPS
14 c
15 c ==================================================================
16 c SUBROUTINE dgscale
17 c ==================================================================
18
19 implicit none
20
21 #include <blas1.h>
22
23 integer nn
24 double precision gold(nn), xdiff(nn), diag(nn)
25
26 integer i
27 double precision r1, rmin, den
28
29
30 c-----------------------------------------
31 c read diagonal
32 c-----------------------------------------
33 call dostore( nn, diag, .false., 3 )
34
35 r1 = 0.
36 do i = 1, nn
37 r1 = r1 + gold(i)*gold(i)*diag(i)
38 end do
39 r1 = 1.0 / r1
40
41 call SSCAL( nn, r1, diag, 1 )
42
43 c-----------------------------------------
44 c update the diagonal
45 c (gg is used as an auxiliary vector)
46 c-----------------------------------------
47
48 den = 0.0
49
50 do i = 1, nn
51 cph(
52 if (diag(i).LE.0) then
53 print *, 'pathei-lsopt: in dgscale; diag = 0 for i=', i
54 diag(i) = rmin
55 end if
56 cph)
57 den = den + xdiff(i)*xdiff(i) / diag(i)
58 end do
59
60 do i = 1, nn
61 diag(i) = 1./
62 $ (1./diag(i)+gold(i)**2-(xdiff(i)/diag(i))**2/den)
63 if (diag(i).le.0.) then
64 diag(i) = rmin
65 endif
66 end do
67
68 c-----------------------------------------
69 c write diagonal
70 c-----------------------------------------
71 call dostore( nn, diag, .true., 3 )
72
73 return
74 end

  ViewVC Help
Powered by ViewVC 1.1.22