ROMS
Loading...
Searching...
No Matches
propagator_so_semi.h
Go to the documentation of this file.
1 MODULE propagator_mod
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! Stochastic Optimals, Seminorm Estimation: !
11! !
12! This routine is used during the computation of the eigenvectors of !
13! the stochastic optimals operator with respect the seminorm of the !
14! chosen functional. !
15! !
16! Reference: !
17! !
18! Moore, A.M. et al., 2004: A comprehensive ocean prediction and !
19! analysis system based on the tangent linear and adjoint of a !
20! regional ocean model, Ocean Modelling, 7, 227-258. !
21! !
22!=======================================================================
23!
24 USE mod_kinds
25!
26 implicit none
27!
28 PRIVATE
29 PUBLIC :: propagator_so_semi
30!
31 CONTAINS
32!
33!***********************************************************************
34 SUBROUTINE propagator_so_semi (RunInterval, state, ad_state)
35!***********************************************************************
36!
37 USE mod_param
38 USE mod_parallel
39 USE mod_iounits
40 USE mod_ocean
41 USE mod_scalars
42 USE mod_stepping
43!
44#ifdef SO_SEMI_WHITE
45 USE packing_mod, ONLY : so_semi_white
46#else
47 USE packing_mod, ONLY : so_semi_red
48#endif
49 USE strings_mod, ONLY : founderror
50!
51! Imported variable declarations.
52!
53 real(dp), intent(in) :: runinterval
54!
55 TYPE (t_gst), intent(in) :: state(ngrids)
56 TYPE (t_gst), intent(inout) :: ad_state(ngrids)
57!
58! Local variable declarations.
59!
60 logical, save :: firstpass = .true.
61!
62 integer :: ng, tile
63!
64 character (len=*), parameter :: myfile = &
65 & __FILE__
66!
67!=======================================================================
68! Backward integration of adjoint model forced with the seminorm of
69! the chosen functional. The adjoint model is run only only once in
70! the first iteration.
71!=======================================================================
72!
73 nrun=nrun+1
74 DO ng=1,ngrids
75 sorec(ng)=0
76 END DO
77!
78 first_pass : IF (firstpass) THEN
79 firstpass=.false.
80!
81! Initialize the adjoint model always from rest.
82!
83 DO ng=1,ngrids
84 CALL ad_initial (ng)
85 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
86 END DO
87!
88! Activate adjoint output.
89!
90 DO ng=1,ngrids
91 ldefadj(ng)=.true.
92 lwrtadj(ng)=.true.
93 lcycleadj(ng)=.false.
94 END DO
95!
96! Time-step adjoint model forced with chosen functional at initial
97! time only.
98!
99 DO ng=1,ngrids
100 IF (master) THEN
101 WRITE (stdout,10) 'AD', ng, ntstart(ng), ntend(ng)
102 END IF
103 dstrs(ng)=tdays(ng)
104 dends(ng)=dstrs(ng)
105 END DO
106
107#ifdef SOLVE3D
108 CALL ad_main3d (runinterval)
109#else
110 CALL ad_main2d (runinterval)
111#endif
112 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
113
114 END IF first_pass
115!
116!-----------------------------------------------------------------------
117! Compute new packed adjoint state vector containing surface forcing
118! variables.
119!-----------------------------------------------------------------------
120!
121 DO ng=1,ngrids
122 DO tile=first_tile(ng),last_tile(ng),+1
123# ifdef SO_SEMI_WHITE
124 CALL so_semi_white (ng, tile, nstr(ng), nend(ng), &
125 & state(ng)%vector, &
126 & ad_state(ng)%vector)
127# else
128 CALL so_semi_red (ng, tile, nstr(ng), nend(ng), &
129 & state(ng)%vector, &
130 & ad_state(ng)%vector)
131# endif
132 END DO
133 END DO
134!
135!-----------------------------------------------------------------------
136! Report iteration and trace or stochastic optimals matrix.
137!-----------------------------------------------------------------------
138!
139 IF (master) THEN
140 DO ng=1,ngrids
141 WRITE (stdout,20) ' PROPAGATOR - Grid: ', ng, &
142 & ', Iteration: ', nrun, &
143 & ', number converged RITZ values: ', &
144 & nconv(ng), 'TRnorm = ', trnorm(ng)
145 END DO
146 END IF
147!
148 10 FORMAT (/,1x,a,1x,'ROMS: started time-stepping:', &
149 & ' (Grid: ',i2.2,' TimeSteps: ',i8.8,' - ',i8.8,')')
150 20 FORMAT (/,a,i2.2,a,i3.3,a,i3.3,/,42x,a,1p,e15.8)
151!
152 RETURN
153 END SUBROUTINE propagator_so_semi
154
155 END MODULE propagator_mod
subroutine ad_initial(ng)
Definition ad_initial.F:4
subroutine ad_main2d
Definition ad_main2d.F:586
subroutine ad_main3d(runinterval)
Definition ad_main3d.F:4
integer stdout
integer, parameter dp
Definition mod_kinds.F:25
integer, dimension(:), allocatable first_tile
logical master
integer, dimension(:), allocatable last_tile
integer, dimension(:), allocatable nstr
Definition mod_param.F:646
integer ngrids
Definition mod_param.F:113
integer, dimension(:), allocatable nend
Definition mod_param.F:647
integer, dimension(:), allocatable sorec
integer, dimension(:), allocatable nconv
real(dp), dimension(:), allocatable tdays
real(r8), dimension(:), allocatable dends
logical, dimension(:), allocatable ldefadj
logical, dimension(:), allocatable lcycleadj
logical, dimension(:), allocatable lwrtadj
integer, dimension(:), allocatable ntend
integer exit_flag
real(r8), dimension(:), allocatable dstrs
integer, dimension(:), allocatable ntstart
integer nrun
real(r8), dimension(:), allocatable trnorm
integer noerror
subroutine, public propagator_so_semi(runinterval, state, ad_state)
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine so_semi_white(ng, tile, mstr, mend, state, ad_state)
Definition packing.F:8788
subroutine so_semi_red(ng, tile, mstr, mend, state, ad_state)
Definition packing.F:8981