SCALE-RM
Functions/Subroutines | Variables
scale_spnudge Module Reference

Functions/Subroutines

subroutine, public spnudge_setup (KA, KS, KE, IA, IS, IE, JA, JS, JE)
 
subroutine, public spnudge_finalize
 

Variables

logical, public spnudge_uv = .false.
 
logical, public spnudge_uv_divfree = .false.
 
integer, public spnudge_uv_lm = 3
 
integer, public spnudge_uv_mm = 3
 
logical, public spnudge_pt = .false.
 
integer, public spnudge_pt_lm = 3
 
integer, public spnudge_pt_mm = 3
 
logical, public spnudge_qv = .false.
 
integer, public spnudge_qv_lm = 3
 
integer, public spnudge_qv_mm = 3
 
real(rp), dimension(:,:,:), allocatable, public spnudge_u_alpha
 
real(rp), dimension(:,:,:), allocatable, public spnudge_v_alpha
 
real(rp), dimension(:,:,:), allocatable, public spnudge_pt_alpha
 
real(rp), dimension(:,:,:), allocatable, public spnudge_qv_alpha
 

Detailed Description

NAMELIST
  • PARAM_SPNUDGE
    nametypedefault valuecomment
    SPNUDGE_UV logical .false.
    SPNUDGE_UV_DIVFREE logical .false.
    SPNUDGE_UV_LM integer 3
    SPNUDGE_UV_MM integer 3
    SPNUDGE_UV_TAU real(RP) 0.0_RP
    SPNUDGE_PT logical .false.
    SPNUDGE_PT_LM integer 3
    SPNUDGE_PT_MM integer 3
    SPNUDGE_PT_TAU real(RP) 0.0_RP
    SPNUDGE_QV logical .false.
    SPNUDGE_QV_LM integer 3
    SPNUDGE_QV_MM integer 3
    SPNUDGE_QV_TAU real(RP) 0.0_RP
    SPNUDGE_LEVEL1 real(RP) 0.0_RP > alpha = 0 for z < SPNUDGE_level1
    SPNUDGE_LEVEL2 real(RP) 0.0_RP
    SPNUDGE_LEVEL3 real(RP) 1E10_RP > alpha = alpha for SPNUDGE_level2 <= z < SPNUDGE_level3
    SPNUDGE_LEVEL4 real(RP) 1E10_RP > alpha = 0 for z => SPNUDGE_level4

History Output
No history output

Function/Subroutine Documentation

◆ spnudge_setup()

subroutine, public scale_spnudge::spnudge_setup ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE 
)

Definition at line 51 of file scale_spnudge.F90.

51  use scale_dft, only: &
52  dft_setup
54  real_cz => atmos_grid_cartesc_real_cz, &
55  real_czuy => atmos_grid_cartesc_real_czuy, &
56  real_czxv => atmos_grid_cartesc_real_czxv
57  implicit none
58  integer, intent(in) :: KA, KS, KE
59  integer, intent(in) :: IA, IS, IE
60  integer, intent(in) :: JA, JS, JE
61 
62  real(RP) :: SPNUDGE_uv_tau = 0.0_rp
63  real(RP) :: SPNUDGE_pt_tau = 0.0_rp
64  real(RP) :: SPNUDGE_qv_tau = 0.0_rp
65 
66  real(RP) :: SPNUDGE_level1 = 0.0_rp
67  real(RP) :: SPNUDGE_level2 = 0.0_rp
68  real(RP) :: SPNUDGE_level3 = 1e10_rp
69  real(RP) :: SPNUDGE_level4 = 1e10_rp
70 
71  namelist /param_spnudge/ &
72  spnudge_uv, &
74  spnudge_uv_lm, &
75  spnudge_uv_mm, &
76  spnudge_uv_tau, &
77  spnudge_pt, &
78  spnudge_pt_lm, &
79  spnudge_pt_mm, &
80  spnudge_pt_tau, &
81  spnudge_qv, &
82  spnudge_qv_lm, &
83  spnudge_qv_mm, &
84  spnudge_qv_tau, &
85  spnudge_level1, &
86  spnudge_level2, &
87  spnudge_level3, &
88  spnudge_level4
89 
90  real(RP) :: uv_alpha, pt_alpha, qv_alpha
91 
92 
93  integer :: k, i, j
94  integer :: ierr
95  !---------------------------------------------------------------------------
96 
97  log_newline
98  log_info("SPNUDGE_setup",*) 'Setup'
99 
100  !--- read namelist
101  rewind(io_fid_conf)
102  read(io_fid_conf,nml=param_spnudge,iostat=ierr)
103  if( ierr < 0 ) then !--- missing
104  log_info("SPNUDGE_setup",*) 'Not found namelist. Default used.'
105  elseif( ierr > 0 ) then !--- fatal error
106  log_error("SPNUDGE_setup",*) 'Not appropriate names in namelist PARAM_SPNUDGE. Check!'
107  call prc_abort
108  endif
109  log_nml(param_spnudge)
110 
111  if ( spnudge_level1 > spnudge_level2 ) then
112  log_error("SPNUDGE_setup",*) 'SPNUDGE_level1 must be lowere or equal to SPNUDGE_level2'
113  call prc_abort
114  end if
115  if ( spnudge_level2 > spnudge_level3 ) then
116  log_error("SPNUDGE_setup",*) 'SPNUDGE_level2 must be lowere or equal to SPNUDGE_level3'
117  call prc_abort
118  end if
119  if ( spnudge_level3 > spnudge_level4 ) then
120  log_error("SPNUDGE_setup",*) 'SPNUDGE_level3 must be lowere or equal to SPNUDGE_level4'
121  call prc_abort
122  end if
123 
124  if ( spnudge_uv .or. spnudge_pt .or. spnudge_qv ) then
125  log_warn("SPNUDGE_setup",*) 'Spectrul nudging is still experimental'
126  end if
127 
128  allocate( spnudge_u_alpha(ka,ia,ja) )
129  allocate( spnudge_v_alpha(ka,ia,ja) )
130  allocate( spnudge_pt_alpha(ka,ia,ja) )
131  allocate( spnudge_qv_alpha(ka,ia,ja) )
132 
133  if ( spnudge_uv .and. spnudge_uv_tau > 0.0_rp ) then
134  uv_alpha = 1.0_rp / spnudge_uv_tau
135 
136  !$omp parallel do
137  do j = js, je
138  do i = is, ie
139  do k = ks, ke
140  if ( real_czuy(k,i,j) < spnudge_level1 ) then
141  spnudge_u_alpha(k,i,j) = 0.0_rp
142  elseif ( real_czuy(k,i,j) < spnudge_level2 ) then
143  spnudge_u_alpha(k,i,j) = ( real_czuy(k,i,j) - spnudge_level1 ) / ( spnudge_level2 - spnudge_level1 ) * uv_alpha
144  elseif ( real_czuy(k,i,j) < spnudge_level3 ) then
145  spnudge_u_alpha(k,i,j) = uv_alpha
146  elseif ( real_czuy(k,i,j) < spnudge_level4 ) then
147  spnudge_u_alpha(k,i,j) = ( real_czuy(k,i,j) - spnudge_level3 ) / ( spnudge_level4 - spnudge_level3 ) * uv_alpha
148  else
149  spnudge_u_alpha(k,i,j) = 0.0_rp
150  endif
151  enddo
152  enddo
153  enddo
154 
155  !$omp parallel do
156  do j = js, je
157  do i = is, ie
158  do k = ks, ke
159  if ( real_czxv(k,i,j) < spnudge_level1 ) then
160  spnudge_v_alpha(k,i,j) = 0.0_rp
161  elseif ( real_czxv(k,i,j) < spnudge_level2 ) then
162  spnudge_v_alpha(k,i,j) = ( real_czxv(k,i,j) - spnudge_level1 ) / ( spnudge_level2 - spnudge_level1 ) * uv_alpha
163  elseif ( real_czxv(k,i,j) < spnudge_level3 ) then
164  spnudge_v_alpha(k,i,j) = uv_alpha
165  elseif ( real_czxv(k,i,j) < spnudge_level4 ) then
166  spnudge_v_alpha(k,i,j) = ( real_czxv(k,i,j) - spnudge_level3 ) / ( spnudge_level4 - spnudge_level3 ) * uv_alpha
167  else
168  spnudge_v_alpha(k,i,j) = 0.0_rp
169  endif
170  enddo
171  enddo
172  enddo
173 
174  else
175 
176  !$omp parallel do
177  do j = js, je
178  do i = is, ie
179  do k = ks, ke
180  spnudge_u_alpha(k,i,j) = 0.0_rp
181  spnudge_v_alpha(k,i,j) = 0.0_rp
182  end do
183  end do
184  end do
185 
186  endif
187 
188  if ( spnudge_pt .and. spnudge_pt_tau > 0.0_rp ) then
189  pt_alpha = 1.0_rp / spnudge_pt_tau
190 
191  !$omp parallel do
192  do j = js, je
193  do i = is, ie
194  do k = ks, ke
195  if ( real_cz(k,i,j) < spnudge_level1 ) then
196  spnudge_pt_alpha(k,i,j) = 0.0_rp
197  elseif ( real_cz(k,i,j) < spnudge_level2 ) then
198  spnudge_pt_alpha(k,i,j) = ( real_cz(k,i,j) - spnudge_level1 ) / ( spnudge_level2 - spnudge_level1 ) * pt_alpha
199  elseif ( real_cz(k,i,j) < spnudge_level3 ) then
200  spnudge_pt_alpha(k,i,j) = pt_alpha
201  elseif ( real_cz(k,i,j) < spnudge_level4 ) then
202  spnudge_pt_alpha(k,i,j) = ( real_cz(k,i,j) - spnudge_level3 ) / ( spnudge_level4 - spnudge_level3 ) * pt_alpha
203  else
204  spnudge_pt_alpha(k,i,j) = 0.0_rp
205  endif
206  enddo
207  enddo
208  enddo
209 
210  else
211 
212  !$omp parallel do
213  do j = js, je
214  do i = is, ie
215  do k = ks, ke
216  spnudge_pt_alpha(k,i,j) = 0.0_rp
217  end do
218  end do
219  end do
220 
221  endif
222 
223  if ( spnudge_qv .and. spnudge_qv_tau > 0.0_rp ) then
224  qv_alpha = 1.0_rp / spnudge_qv_tau
225 
226  !$omp parallel do
227  do j = js, je
228  do i = is, ie
229  do k = ks, ke
230  if ( real_cz(k,i,j) < spnudge_level1 ) then
231  spnudge_qv_alpha(k,i,j) = 0.0_rp
232  elseif ( real_cz(k,i,j) < spnudge_level2 ) then
233  spnudge_qv_alpha(k,i,j) = ( real_cz(k,i,j) - spnudge_level1 ) / ( spnudge_level2 - spnudge_level1 ) * qv_alpha
234  elseif ( real_cz(k,i,j) < spnudge_level3 ) then
235  spnudge_qv_alpha(k,i,j) = qv_alpha
236  elseif ( real_cz(k,i,j) < spnudge_level4 ) then
237  spnudge_qv_alpha(k,i,j) = ( real_cz(k,i,j) - spnudge_level3 ) / ( spnudge_level4 - spnudge_level3 ) * qv_alpha
238  else
239  spnudge_qv_alpha(k,i,j) = 0.0_rp
240  endif
241  enddo
242  enddo
243  enddo
244 
245  else
246 
247  !$omp parallel do
248  do j = js, je
249  do i = is, ie
250  do k = ks, ke
251  spnudge_qv_alpha(k,i,j) = 0.0_rp
252  end do
253  end do
254  end do
255 
256  end if
257 
258  !$acc enter data copyin(SPNUDGE_u_alpha, SPNUDGE_v_alpha, SPNUDGE_pt_alpha, SPNUDGE_qv_alpha)
259 
260  call dft_setup( ka, ks, ke, ia, is, ie, ja, js, je, &
263 
264  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_czuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_czxv, scale_dft::dft_setup(), scale_io::io_fid_conf, scale_prc::prc_abort(), spnudge_pt, spnudge_pt_alpha, spnudge_pt_lm, spnudge_pt_mm, spnudge_qv, spnudge_qv_alpha, spnudge_qv_lm, spnudge_qv_mm, spnudge_u_alpha, spnudge_uv, spnudge_uv_divfree, spnudge_uv_lm, spnudge_uv_mm, and spnudge_v_alpha.

Referenced by scale_atmos_dyn::atmos_dyn_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ spnudge_finalize()

subroutine, public scale_spnudge::spnudge_finalize

Definition at line 268 of file scale_spnudge.F90.

268  use scale_dft, only: &
270 
271  call dft_finalize
272 
273  !$acc exit data delete(SPNUDGE_u_alpha, SPNUDGE_v_alpha, SPNUDGE_pt_alpha, SPNUDGE_qv_alpha)
274  deallocate( spnudge_u_alpha )
275  deallocate( spnudge_v_alpha )
276  deallocate( spnudge_pt_alpha )
277  deallocate( spnudge_qv_alpha )
278 
279  return

References scale_dft::dft_finalize(), spnudge_pt_alpha, spnudge_qv_alpha, spnudge_u_alpha, and spnudge_v_alpha.

Referenced by scale_atmos_dyn::atmos_dyn_finalize().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ spnudge_uv

logical, public scale_spnudge::spnudge_uv = .false.

Definition at line 25 of file scale_spnudge.F90.

25  logical, public :: SPNUDGE_uv = .false.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_uv_divfree

logical, public scale_spnudge::spnudge_uv_divfree = .false.

Definition at line 26 of file scale_spnudge.F90.

26  logical, public :: SPNUDGE_uv_divfree = .false.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_uv_lm

integer, public scale_spnudge::spnudge_uv_lm = 3

Definition at line 27 of file scale_spnudge.F90.

27  integer, public :: SPNUDGE_uv_lm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_uv_mm

integer, public scale_spnudge::spnudge_uv_mm = 3

Definition at line 28 of file scale_spnudge.F90.

28  integer, public :: SPNUDGE_uv_mm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_pt

logical, public scale_spnudge::spnudge_pt = .false.

Definition at line 30 of file scale_spnudge.F90.

30  logical, public :: SPNUDGE_pt = .false.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_pt_lm

integer, public scale_spnudge::spnudge_pt_lm = 3

Definition at line 31 of file scale_spnudge.F90.

31  integer, public :: SPNUDGE_pt_lm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_pt_mm

integer, public scale_spnudge::spnudge_pt_mm = 3

Definition at line 32 of file scale_spnudge.F90.

32  integer, public :: SPNUDGE_pt_mm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_qv

logical, public scale_spnudge::spnudge_qv = .false.

Definition at line 34 of file scale_spnudge.F90.

34  logical, public :: SPNUDGE_qv = .false.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_qv_lm

integer, public scale_spnudge::spnudge_qv_lm = 3

Definition at line 35 of file scale_spnudge.F90.

35  integer, public :: SPNUDGE_qv_lm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_qv_mm

integer, public scale_spnudge::spnudge_qv_mm = 3

Definition at line 36 of file scale_spnudge.F90.

36  integer, public :: SPNUDGE_qv_mm = 3

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), and spnudge_setup().

◆ spnudge_u_alpha

real(rp), dimension(:,:,:), allocatable, public scale_spnudge::spnudge_u_alpha

Definition at line 38 of file scale_spnudge.F90.

38  real(RP), allocatable, public :: SPNUDGE_u_alpha(:,:,:)

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), spnudge_finalize(), and spnudge_setup().

◆ spnudge_v_alpha

real(rp), dimension(:,:,:), allocatable, public scale_spnudge::spnudge_v_alpha

Definition at line 39 of file scale_spnudge.F90.

39  real(RP), allocatable, public :: SPNUDGE_v_alpha(:,:,:)

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), spnudge_finalize(), and spnudge_setup().

◆ spnudge_pt_alpha

real(rp), dimension(:,:,:), allocatable, public scale_spnudge::spnudge_pt_alpha

Definition at line 40 of file scale_spnudge.F90.

40  real(RP), allocatable, public :: SPNUDGE_pt_alpha(:,:,:)

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), spnudge_finalize(), and spnudge_setup().

◆ spnudge_qv_alpha

real(rp), dimension(:,:,:), allocatable, public scale_spnudge::spnudge_qv_alpha

Definition at line 41 of file scale_spnudge.F90.

41  real(RP), allocatable, public :: SPNUDGE_qv_alpha(:,:,:)

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), spnudge_finalize(), and spnudge_setup().

scale_dft::dft_finalize
subroutine, public dft_finalize
Definition: scale_dft.F90:102
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_spnudge::spnudge_uv_lm
integer, public spnudge_uv_lm
Definition: scale_spnudge.F90:27
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
scale_spnudge::spnudge_pt_lm
integer, public spnudge_pt_lm
Definition: scale_spnudge.F90:31
scale_spnudge::spnudge_pt
logical, public spnudge_pt
Definition: scale_spnudge.F90:30
scale_spnudge::spnudge_pt_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_pt_alpha
Definition: scale_spnudge.F90:40
scale_spnudge::spnudge_uv_divfree
logical, public spnudge_uv_divfree
Definition: scale_spnudge.F90:26
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_spnudge::spnudge_qv_lm
integer, public spnudge_qv_lm
Definition: scale_spnudge.F90:35
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_czxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_czxv
geopotential height [m] (zxv)
Definition: scale_atmos_grid_cartesC_real.F90:41
scale_spnudge::spnudge_qv
logical, public spnudge_qv
Definition: scale_spnudge.F90:34
scale_dft
Definition: scale_dft.F90:3
scale_spnudge::spnudge_pt_mm
integer, public spnudge_pt_mm
Definition: scale_spnudge.F90:32
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_spnudge::spnudge_u_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_u_alpha
Definition: scale_spnudge.F90:38
scale_dft::dft_setup
subroutine, public dft_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE, LM, MM)
Definition: scale_dft.F90:37
scale_spnudge::spnudge_qv_mm
integer, public spnudge_qv_mm
Definition: scale_spnudge.F90:36
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_spnudge::spnudge_uv_mm
integer, public spnudge_uv_mm
Definition: scale_spnudge.F90:28
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_spnudge::spnudge_qv_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_qv_alpha
Definition: scale_spnudge.F90:41
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_spnudge::spnudge_uv
logical, public spnudge_uv
Definition: scale_spnudge.F90:25
scale_spnudge::spnudge_v_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_v_alpha
Definition: scale_spnudge.F90:39
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_czuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_czuy
geopotential height [m] (zuy)
Definition: scale_atmos_grid_cartesC_real.F90:40