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 branches/DEV_1879_mpp_rep/NEMO/OPA_SRC/SOL – NEMO

source: branches/DEV_1879_mpp_rep/NEMO/OPA_SRC/SOL/solpcg.F90 @ 1920

Last change on this file since 1920 was 1920, checked in by rblod, 14 years ago

Add modifications for mpp reproducibility, see ticket #677

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 8.6 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   USE oce             ! ocean dynamics and tracers variables
11   USE dom_oce         ! ocean space and time domain variables
12   USE sol_oce         ! ocean solver variables
13   USE lib_mpp         ! distributed memory computing
14   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
15   USE in_out_manager  ! I/O manager
16   USE lib_fortran
17
18   IMPLICIT NONE
19   PRIVATE
20
21   PUBLIC   sol_pcg    !
22
23   !! * Substitutions
24#  include "vectopt_loop_substitute.h90"
25   !!----------------------------------------------------------------------
26   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)
27   !! $Id$
28   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
29   !!----------------------------------------------------------------------
30CONTAINS
31
32   SUBROUTINE sol_pcg( kindic )
33      !!----------------------------------------------------------------------
34      !!                  ***  ROUTINE sol_pcg  ***
35      !!                   
36      !! ** Purpose :   Solve the ellipic equation for the transport
37      !!      divergence system  using a diagonal preconditionned
38      !!      conjugate gradient method.
39      !!
40      !! ** Method  :   Diagonal preconditionned conjugate gradient method.
41      !!      the algorithm is multitasked. (case of 5 points matrix)
42      !!      define              pa  = q^-1 * a
43      !!                        pgcb  = q^-1 * gcb
44      !!                 < . ; . >_q  = ( . )^t q ( . )
45      !!      where q is the preconditioning matrix = diagonal matrix of the
46      !!                                              diagonal elements of a
47      !!      Initialization  :
48      !!         x(o) = gcx
49      !!         r(o) = d(o) = pgcb - pa.x(o)
50      !!         rr(o)= < r(o) , r(o) >_q
51      !!      Iteration 1     :
52      !!         standard PCG algorithm
53      !!      Iteration n > 1 :
54      !!         s(n)   = pa.r(n)
55      !!         gam(n) = < r(n) , r(n) >_q
56      !!         del(n) = < r(n) , s(n) >_q
57      !!         bet(n) = gam(n) / gam(n-1)
58      !!         d(n)   = r(n) + bet(n) d(n-1)
59      !!         z(n)   = s(n) + bet(n) z(n-1)
60      !!         sig(n) = del(n) - bet(n)*bet(n)*sig(n-1)
61      !!         alp(n) = gam(n) / sig(n)
62      !!         x(n+1) = x(n) + alp(n) d(n)
63      !!         r(n+1) = r(n) - alp(n) z(n)
64      !!      Convergence test :
65      !!         rr(n+1) / < gcb , gcb >_q   =< epsr
66      !!
67      !! ** Action : - niter  : solver number of iteration done
68      !!             - res    : solver residu reached
69      !!             - gcx()  : solution of the elliptic system
70      !!
71      !! References :
72      !!      Madec et al. 1988, Ocean Modelling, issue 78, 1-6.
73      !!      D Azevedo et al. 1993, Computer Science Technical Report, Tennessee U.
74      !!
75      !! History :
76      !!        !  90-10  (G. Madec)  Original code
77      !!        !  91-11  (G. Madec)
78      !!        !  93-04  (M. Guyon)  loops and suppress pointers
79      !!        !  95-09  (M. Imbard, J. Escobar)  mpp exchange
80      !!        !  96-05  (G. Madec)  merge sor and pcg formulations
81      !!        !  96-11  (A. Weaver)  correction to preconditioning
82      !!   8.5  !  02-08  (G. Madec)  F90: Free form
83      !!        !  08-01  (R. Benshila) mpp optimization
84      !!----------------------------------------------------------------------
85      INTEGER, INTENT( inout ) ::   kindic   ! solver indicator, < 0 if the conver-
86      !                                      ! gence is not reached: the model is
87      !                                      ! stopped in step
88      !                                      ! set to zero before the call of solpcg
89      !!
90      INTEGER ::   ji, jj, jn                ! dummy loop indices
91      REAL(wp) ::  zgcad                     ! temporary scalars
92      REAL(wp), DIMENSION(2) :: zsum
93      REAL(wp), DIMENSION(jpi,jpj) :: zgcr
94      !!----------------------------------------------------------------------
95
96      ! Initialization of the algorithm with standard PCG
97      ! -------------------------------------------------
98      zgcr = 0.e0 
99      gcr  = 0.e0 
100
101      CALL lbc_lnk( gcx, c_solver_pt, 1. )   ! lateral boundary condition
102
103      ! gcr   = gcb-a.gcx
104      ! gcdes = gcr
105      DO jj = 2, jpjm1
106         DO ji = fs_2, fs_jpim1   ! vector opt.
107            zgcad = bmask(ji,jj) * ( gcb(ji,jj  ) -                gcx(ji  ,jj  )   &
108               &                                  - gcp(ji,jj,1) * gcx(ji  ,jj-1)   &
109               &                                  - gcp(ji,jj,2) * gcx(ji-1,jj  )   &
110               &                                  - gcp(ji,jj,3) * gcx(ji+1,jj  )   &
111               &                                  - gcp(ji,jj,4) * gcx(ji  ,jj+1)   )
112            gcr  (ji,jj) = zgcad
113            gcdes(ji,jj) = zgcad
114         END DO
115      END DO
116
117      ! rnorme = (gcr,gcr)
118      rnorme = glob_sum(  gcr(:,:) * gcdmat(:,:) * gcr(:,:)  )
119
120      CALL lbc_lnk( gcdes, c_solver_pt, 1. )   ! lateral boundary condition
121
122      ! gccd = matrix . gcdes
123      DO jj = 2, jpjm1
124         DO ji = fs_2, fs_jpim1   ! vector opt.
125            gccd(ji,jj) = bmask(ji,jj)*( gcdes(ji,jj)   &
126               &        +gcp(ji,jj,1)*gcdes(ji,jj-1)+gcp(ji,jj,2)*gcdes(ji-1,jj)   &
127               &        +gcp(ji,jj,4)*gcdes(ji,jj+1)+gcp(ji,jj,3)*gcdes(ji+1,jj)   )
128         END DO
129      END DO 
130
131      ! alph = (gcr,gcr)/(gcdes,gccd)
132      radd = glob_sum(  gcdes(:,:) * gcdmat(:,:) * gccd(:,:)  )
133      alph = rnorme /radd
134
135      ! gcx = gcx + alph * gcdes
136      ! gcr = gcr - alph * gccd
137      DO jj = 2, jpjm1
138         DO ji = fs_2, fs_jpim1   ! vector opt.
139            gcx(ji,jj) = bmask(ji,jj) * ( gcx(ji,jj) + alph * gcdes(ji,jj) )
140            gcr(ji,jj) = bmask(ji,jj) * ( gcr(ji,jj) - alph * gccd (ji,jj) )
141         END DO
142      END DO
143
144      ! Algorithm wtih Eijkhout rearrangement
145      ! -------------------------------------
146       
147      !                                                !================
148      DO jn = 1, nn_nmax                               ! Iterative loop
149         !                                             !================
150
151         CALL lbc_lnk( gcr, c_solver_pt, 1. )   ! lateral boundary condition
152
153         ! zgcr = matrix . gcr
154         DO jj = 2, jpjm1
155            DO ji = fs_2, fs_jpim1   ! vector opt.
156               zgcr(ji,jj) = bmask(ji,jj)*( gcr(ji,jj)   &
157                  &        +gcp(ji,jj,1)*gcr(ji,jj-1)+gcp(ji,jj,2)*gcr(ji-1,jj)   &
158                  &        +gcp(ji,jj,4)*gcr(ji,jj+1)+gcp(ji,jj,3)*gcr(ji+1,jj)   )
159            END DO
160         END DO
161 
162         ! rnorme = (gcr,gcr)
163         rr = rnorme
164
165         ! zgcad = (zgcr,gcr)
166         zsum(1) = glob_sum(gcr(:,:) * gcdmat(:,:) * gcr(:,:))
167         zsum(2) = glob_sum(gcr(:,:) * gcdmat(:,:) * zgcr(:,:) * bmask(:,:))
168
169         !!RB we should gather the 2 glob_sum
170         rnorme = zsum(1) 
171         zgcad  = zsum(2)
172         ! test of convergence
173         IF( rnorme < epsr .OR. jn == nn_nmax ) THEN
174            res = SQRT( rnorme )
175            niter = jn
176            ncut = 999
177         ENDIF
178
179         ! beta = (rk+1,rk+1)/(rk,rk)
180         beta = rnorme / rr
181         radd = zgcad - beta*beta*radd
182         alph = rnorme / radd
183
184         ! gcx = gcx + alph * gcdes
185         ! gcr = gcr - alph * gccd
186         DO jj = 2, jpjm1
187            DO ji = fs_2, fs_jpim1   ! vector opt.
188               gcdes(ji,jj) = gcr (ji,jj) + beta * gcdes(ji,jj) 
189               gccd (ji,jj) = zgcr(ji,jj) + beta * gccd (ji,jj) 
190               gcx  (ji,jj) = gcx (ji,jj) + alph * gcdes(ji,jj) 
191               gcr  (ji,jj) = gcr (ji,jj) - alph * gccd (ji,jj) 
192            END DO
193         END DO
194       
195         ! indicator of non-convergence or explosion
196         IF( jn == nn_nmax .OR. SQRT(epsr)/eps > 1.e+20 ) kindic = -2
197         IF( ncut == 999 ) GOTO 999
198
199         !                                             !================
200      END DO                                           !    End Loop
201      !                                                !================
202     
203999   CONTINUE
204     
205     
206      ! Output in gcx with lateral b.c. applied
207      ! ---------------------------------------
208     
209      CALL lbc_lnk( gcx, c_solver_pt, 1. )
210     
211   END SUBROUTINE sol_pcg
212
213   !!=====================================================================
214END MODULE solpcg
Note: See TracBrowser for help on using the repository browser.