New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
solpcg.F90 in trunk/NEMO/OPA_SRC/SOL – NEMO

source: trunk/NEMO/OPA_SRC/SOL/solpcg.F90 @ 16

Last change on this file since 16 was 16, checked in by opalod, 20 years ago

CT : UPDATE001 : First major NEMO update

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 7.8 KB
Line 
1MODULE solpcg
2   !!======================================================================
3   !!                     ***  MODULE  solfet
4   !! Ocean solver :  preconditionned conjugate gradient solver
5   !!=====================================================================
6
7   !!----------------------------------------------------------------------
8   !!   sol_pcg    : preconditionned conjugate gradient solver
9   !!----------------------------------------------------------------------
10   !! * Modules used
11   USE oce             ! ocean dynamics and tracers variables
12   USE dom_oce         ! ocean space and time domain variables
13   USE sol_oce         ! ocean solver variables
14   USE lib_mpp         ! distributed memory computing
15   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
16   USE in_out_manager  ! I/O manager
17
18   IMPLICIT NONE
19   PRIVATE
20
21   !! * Routine accessibility
22   PUBLIC sol_pcg              ! ???
23
24   !! * Substitutions
25#  include "vectopt_loop_substitute.h90"
26   !!----------------------------------------------------------------------
27
28CONTAINS
29
30   SUBROUTINE sol_pcg( kindic )
31      !!----------------------------------------------------------------------
32      !!                  ***  ROUTINE sol_pcg  ***
33      !!                   
34      !! ** Purpose :   Solve the ellipic equation for the barotropic stream
35      !!      function system (lk_dynspg_rl=T) or the transport divergence
36      !!      system (lk_dynspg_fsc=T) using a diagonal preconditionned
37      !!      conjugate gradient method.
38      !!      In the former case, the barotropic stream function trend has a
39      !!      zero boundary condition along all coastlines (i.e. continent
40      !!      as well as islands) while in the latter the boundary condition
41      !!      specification is not required.
42      !!
43      !! ** Method  :   Diagonal preconditionned conjugate gradient method.
44      !!      the algorithm is multitasked. (case of 5 points matrix)
45      !!      define              pa  = q^-1 * a
46      !!                        pgcb  = q^-1 * gcb
47      !!                 < . ; . >_q  = ( . )^t q ( . )
48      !!      where q is the preconditioning matrix = diagonal matrix of the
49      !!                                              diagonal elements of a
50      !!      Initialization:
51      !!         x(o) = gcx
52      !!         r(o) = d(o) = pgcb - pa.x(o)
53      !!         rr(o)= < r(o) , r(o) >_q
54      !!      Iteration n   :
55      !!         z(n)   = pa.d(n)
56      !!         alp(n) = rr(n) / < z(n) , d(n) >_q
57      !!         x(n+1) = x(n) + alp(n) d(n)
58      !!         r(n+1) = r(n) - alp(n) z(n)
59      !!         rr(n+1)= < r(n+1) , r(n+1) >_q
60      !!         bet(n) = rr(n+1) / rr(n)
61      !!         r(n+1) = r(n+1) + bet(n+1) d(n)
62      !!      Convergence test :
63      !!         rr(n+1) / < gcb , gcb >_q   =< epsr
64      !!
65      !! ** Action : - niter  : solver number of iteration done
66      !!             - res    : solver residu reached
67      !!             - gcx()  : solution of the elliptic system
68      !!
69      !! References :
70      !!      Madec et al. 1988, Ocean Modelling, issue 78, 1-6.
71      !!
72      !! History :
73      !!        !  90-10  (G. Madec)  Original code
74      !!        !  91-11  (G. Madec)
75      !!        !  93-04  (M. Guyon)  loops and suppress pointers
76      !!        !  95-09  (M. Imbard, J. Escobar)  mpp exchange
77      !!        !  96-05  (G. Madec)  merge sor and pcg formulations
78      !!        !  96-11  (A. Weaver)  correction to preconditioning
79      !!   8.5  !  02-08  (G. Madec)  F90: Free form
80      !!----------------------------------------------------------------------
81      !! * Arguments
82      INTEGER, INTENT( inout ) ::   kindic   ! solver indicator, < 0 if the conver-
83      !                                      ! gence is not reached: the model is
84      !                                      ! stopped in step
85      !                                      ! set to zero before the call of solpcg
86
87      !! * Local declarations
88      INTEGER ::   ji, jj, jn                ! dummy loop indices
89      REAL(wp) ::   zgcad                    ! temporary scalars
90      !!----------------------------------------------------------------------
91
92      !                                                !================
93      DO jn = 1, nmax                                  ! Iterative loop
94         !                                             !================
95
96         IF( jn == 1 ) THEN           ! Initialization of the algorithm
97
98            CALL lbc_lnk( gcx, c_solver_pt, 1. )   ! lateral boundary condition
99
100            ! gcr   = gcb-a.gcx
101            ! gcdes = gsr
102            DO jj = 2, jpjm1
103               DO ji = fs_2, fs_jpim1   ! vector opt.
104                  zgcad = bmask(ji,jj) * ( gcb(ji,jj  ) -                gcx(ji  ,jj  )   &
105                     &                                  - gcp(ji,jj,1) * gcx(ji  ,jj-1)   &
106                     &                                  - gcp(ji,jj,2) * gcx(ji-1,jj  )   &
107                     &                                  - gcp(ji,jj,3) * gcx(ji+1,jj  )   &
108                     &                                  - gcp(ji,jj,4) * gcx(ji  ,jj+1)   )
109                  gcr  (ji,jj) = zgcad
110                  gcdes(ji,jj) = zgcad
111               END DO
112            END DO
113           
114            rnorme = SUM(  gcr(:,:) * gcdmat(:,:) * gcr(:,:)  )
115            IF( lk_mpp )   CALL mpp_sum( rnorme )   ! sum over the global domain
116            rr = rnorme
117
118         ENDIF
119       
120         !                             ! Algorithm
121       
122         CALL lbc_lnk( gcdes, c_solver_pt, 1. )   ! lateral boundary condition
123       
124         ! ... gccd = matrix . gcdes
125         DO jj = 2, jpjm1
126            DO ji = fs_2, fs_jpim1   ! vector opt.
127               gccd(ji,jj) = bmask(ji,jj)*( gcdes(ji,jj)   &
128                  &        +gcp(ji,jj,1)*gcdes(ji,jj-1)+gcp(ji,jj,2)*gcdes(ji-1,jj)   &
129                  &        +gcp(ji,jj,4)*gcdes(ji,jj+1)+gcp(ji,jj,3)*gcdes(ji+1,jj)   )
130            END DO
131         END DO
132 
133         ! alph = (gcr,gcr)/(gcdes,gccd)
134         radd = SUM(  gcdes(:,:) * gcdmat(:,:) * gccd(:,:)  )
135         IF( lk_mpp )   CALL mpp_sum( radd )   ! sum over the global domain
136         alph = rr / radd
137         
138         ! gcx = gcx + alph * gcdes
139         ! gcr = gcr - alph * gccd
140         DO jj = 2, jpjm1
141            DO ji = fs_2, fs_jpim1   ! vector opt.
142               gcx(ji,jj) = bmask(ji,jj) * ( gcx(ji,jj) + alph * gcdes(ji,jj) )
143               gcr(ji,jj) = bmask(ji,jj) * ( gcr(ji,jj) - alph * gccd (ji,jj) )
144            END DO
145         END DO
146       
147         ! rnorme = (gcr,gcr)
148         rnorme = SUM(  gcr(:,:) * gcdmat(:,:) * gcr(:,:)  )
149         IF( lk_mpp )   CALL  mpp_sum( rnorme )   ! sum over the global domain
150       
151         ! test of convergence
152         IF( rnorme < epsr .OR. jn == nmax ) THEN
153            res = SQRT( rnorme )
154            niter = jn
155            ncut = 999
156         ENDIF
157       
158         ! beta = (rk+1,rk+1)/(rk,rk)
159         beta = rnorme / rr
160         rr   = rnorme
161
162         ! indicator of non-convergence or explosion
163         IF( jn == nmax .OR. SQRT(epsr)/eps > 1.e+20 ) kindic = -2
164         IF( ncut == 999 ) GOTO 999
165
166         ! gcdes = gcr + beta * gcdes
167         DO jj = 2, jpjm1
168            DO ji = fs_2, fs_jpim1   ! vector opt.
169               gcdes(ji,jj) = bmask(ji,jj)*( gcr(ji,jj) + beta * gcdes(ji,jj) )
170            END DO
171         END DO
172       
173         !                                             !================
174      END DO                                           !    End Loop
175      !                                                !================
176     
177999   CONTINUE
178     
179     
180      ! Output in gcx with lateral b.c. applied
181      ! ---------------------------------------
182     
183      CALL lbc_lnk( gcx, c_solver_pt, 1. )
184     
185   END SUBROUTINE sol_pcg
186
187   !!=====================================================================
188END MODULE solpcg
Note: See TracBrowser for help on using the repository browser.