SCALE-RM
Data Types | Functions/Subroutines
scale_atmos_dyn_tinteg_rkcommon Module Reference

module scale_atmos_dyn_tinteg_rkcommon More...

Data Types

type  rkinfo
 

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_rkcommon_setup (this, rk_stage_num, rk_register_num, rkcoef_a, rkcoef_b, varname_list, is_type_flux, alloc_rkwork_flag, comm_id_offset)
 
subroutine, public atmos_dyn_tinteg_rkcommon_rkwork_alloc (this)
 
subroutine, public atmos_dyn_tinteg_rkcommon_rkwork_dealloc (this)
 
subroutine, public atmos_dyn_tinteg_rkcommon_comm (this)
 
subroutine, public atmos_dyn_tinteg_rkcommon_comm_wait (this)
 
subroutine, public atmos_dyn_tinteg_rkcommon_nextstage (this, nowstage, io, jo, ko, dt)
 
subroutine, public atmos_dyn_tinteg_rkcommon_updatevar (this, io, jo, ko, vs, ve, dt, var)
 
subroutine, public atmos_dyn_tinteg_rkcommon_updateflux (this, nowstage, io, jo, ko, va_, flux)
 

Detailed Description

module scale_atmos_dyn_tinteg_rkcommon

Description
A modoule providing utility for some runge-kutta schemes
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_rkcommon_setup()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_setup ( type(rkinfo), intent(inout)  this,
integer, intent(in)  rk_stage_num,
integer, intent(in)  rk_register_num,
real(rp), dimension(rk_stage_num,rk_stage_num), intent(in)  rkcoef_a,
real(rp), dimension(rk_stage_num), intent(in)  rkcoef_b,
character(*), dimension(:), intent(in)  varname_list,
logical, intent(in), optional  is_type_flux,
logical, intent(in), optional  alloc_rkwork_flag,
integer, intent(in), optional  comm_id_offset 
)

Definition at line 81 of file scale_atmos_dyn_tinteg_rkcommon.F90.

81 
83  implicit none
84 
85  type(RKInfo), intent(inout) :: this
86  integer, intent(in) :: rk_stage_num
87  integer, intent(in) :: rk_register_num
88  real(RP), intent(in) :: rkcoef_a(rk_stage_num,rk_stage_num)
89  real(RP), intent(in) :: rkcoef_b(rk_stage_num)
90  character(*), intent(in) :: varname_list(:)
91  logical, intent(in), optional :: is_type_flux
92  logical, intent(in), optional :: alloc_rkwork_flag
93  integer, intent(in), optional :: comm_id_offset
94 
95  integer :: reg_id
96  integer :: var_id
97  character(H_MID) :: rktag
98  !------------------------------------------
99 
100  this%var_num = size(varname_list)
101  this%stage_num = rk_stage_num
102  this%register_num = rk_register_num
103  if (present(is_type_flux)) then
104  this%flux_flag = is_type_flux
105  else
106  this%flux_flag = .false.
107  end if
108 
109 
110  allocate( this%rkcoef_a(this%stage_num,this%stage_num) )
111  allocate( this%rkcoef_b(this%stage_num) )
112  this%rkcoef_a(:,:) = rkcoef_a(:,:)
113  this%rkcoef_b(:) = rkcoef_b(:)
114 
115  ! If necessary, allocate working arrays to register the tendencies at some RK stagaes.
116 
117  if ( present(alloc_rkwork_flag) ) then
118  if (alloc_rkwork_flag) call atmos_dyn_tinteg_rkcommon_rkwork_alloc( this )
119  end if
120 
121  ! Allocate an array to store variables at next RK stage
122 
123  allocate( this%buf(ka,ia,ja,this%var_num) )
124 
125  !$omp parallel workshare
126  this%buf(:,:,:,:) = undef
127  !$omp end parallel workshare
128 
129  ! Prepair some information for communication with MPI
130 
131  if ( .not. this%flux_flag ) then
132  allocate( this%comm_ind(this%var_num) )
133  do var_id = 1, this%var_num
134  this%comm_ind(var_id) = var_id
135  end do
136  if (present(comm_id_offset)) &
137  this%comm_ind(:) = comm_id_offset + this%comm_ind(:)
138 
139  do var_id = 1, this%var_num
140  write(rktag,'(a,a)') trim(varname_list(var_id)), 'RK'
141  call comm_vars8_init( trim(rktag), this%buf(:,:,:,var_id), this%comm_ind(var_id) )
142  end do
143  end if
144 
145  return

References atmos_dyn_tinteg_rkcommon_rkwork_alloc(), scale_comm_cartesc::comm_vars8_init(), scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, and scale_atmos_grid_cartesc_index::ka.

Here is the call graph for this function:

◆ atmos_dyn_tinteg_rkcommon_rkwork_alloc()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_rkwork_alloc ( type(rkinfo), intent(inout)  this)

Definition at line 149 of file scale_atmos_dyn_tinteg_rkcommon.F90.

149  implicit none
150 
151  type(RKInfo), intent(inout) :: this
152  !------------------------------------------
153 
154  if (this%register_num > 0 .and. (.not. this%flux_flag) ) then
155 
156  allocate( this%work0(ka,ia,ja,this%var_num) )
157  allocate( this%work(ka,ia,ja,this%var_num,this%register_num) )
158 
159 #ifdef DEBUG
160  !$omp parallel workshare
161  this%work0 (:,:,:,:) = undef
162  this%work (:,:,:,:,:) = undef
163  !$omp end parallel workshare
164 #endif
165  end if
166 
167  return

References scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, and scale_atmos_grid_cartesc_index::ka.

Referenced by atmos_dyn_tinteg_rkcommon_setup().

Here is the caller graph for this function:

◆ atmos_dyn_tinteg_rkcommon_rkwork_dealloc()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_rkwork_dealloc ( type(rkinfo), intent(inout)  this)

Definition at line 171 of file scale_atmos_dyn_tinteg_rkcommon.F90.

171  implicit none
172 
173  type(RKInfo), intent(inout) :: this
174  !------------------------------------------
175 
176  if ( allocated(this%work) ) then
177  deallocate( this%work0 )
178  deallocate( this%work )
179  end if
180 
181  return

◆ atmos_dyn_tinteg_rkcommon_comm()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_comm ( type(rkinfo), intent(inout)  this)

Definition at line 185 of file scale_atmos_dyn_tinteg_rkcommon.F90.

185  use scale_comm_cartesc, only: comm_vars8
186  implicit none
187 
188  type(RKInfo), intent(inout) :: this
189 
190  integer :: var_id
191  !------------------------------------------
192 
193  do var_id = 1, this%var_num
194  call comm_vars8( this%buf(:,:,:,var_id), this%comm_ind(var_id) )
195  end do
196 
197  return

◆ atmos_dyn_tinteg_rkcommon_comm_wait()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_comm_wait ( type(rkinfo), intent(inout)  this)

Definition at line 201 of file scale_atmos_dyn_tinteg_rkcommon.F90.

201  use scale_comm_cartesc, only: comm_wait
202  implicit none
203 
204  type(RKInfo), intent(inout) :: this
205 
206  integer :: var_id
207  !------------------------------------------
208 
209  do var_id = 1, this%var_num
210  call comm_wait( this%buf(:,:,:,var_id), this%comm_ind(var_id), .false. )
211  end do
212 
213  return

◆ atmos_dyn_tinteg_rkcommon_nextstage()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_nextstage ( type(rkinfo), intent(inout)  this,
integer, intent(in)  nowstage,
integer, dimension(this%var_num), intent(in)  io,
integer, dimension(this%var_num), intent(in)  jo,
integer, dimension(this%var_num), intent(in)  ko,
real(rp), intent(in)  dt 
)

Definition at line 217 of file scale_atmos_dyn_tinteg_rkcommon.F90.

217  implicit none
218  type(RKInfo), intent(inout) :: this
219  integer, intent(in) :: nowstage
220  integer, intent(in) :: io(this%var_num)
221  integer, intent(in) :: jo(this%var_num)
222  integer, intent(in) :: ko(this%var_num)
223  real(RP), intent(in) :: dt
224 
225  integer :: i, j, k, iv, rks
226  real(RP) :: var0
227 
228  real(RP) :: a_(this%stage_num)
229  !--------------------------------------
230 
231  a_(:) = dt * this%RKCoef_a(nowstage+1,:)
232 
233  !$omp parallel private(j,i,k,var0,iv,rks)
234 
235  do iv=1, this%var_num
236  !$omp do collapse(2)
237  do j=js-jo(iv), je
238  do i=is-io(iv), ie
239  do k=ks-ko(iv), ke
240  var0 = this%work0(k,i,j,iv)
241  this%work(k,i,j,iv,nowstage) = this%work(k,i,j,iv,nowstage) - var0
242  this%buf(k,i,j,iv) = var0 + a_(nowstage) * this%work(k,i,j,iv,nowstage)
243  end do
244  end do
245  end do
246  !$omp end do
247  end do
248 
249  do rks=1, nowstage-1
250  if ( abs(this%RKCoef_a(nowstage+1,rks)) < eps ) cycle
251  do iv=1, this%var_num
252  !$omp do collapse(2)
253  do j=js-jo(iv), je
254  do i=is-io(iv), ie
255  do k=ks-ko(iv), ke
256  this%buf(k,i,j,iv) = this%buf(k,i,j,iv) + a_(rks) * this%work(k,i,j,iv,rks)
257  end do
258  end do
259  end do
260  !$omp end do
261  end do
262  end do
263 
264  !$omp end parallel
265 
266  return

References scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

◆ atmos_dyn_tinteg_rkcommon_updatevar()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_updatevar ( type(rkinfo), intent(inout)  this,
integer, dimension(this%var_num), intent(in)  io,
integer, dimension(this%var_num), intent(in)  jo,
integer, dimension(this%var_num), intent(in)  ko,
integer, intent(in)  vs,
integer, intent(in)  ve,
real(rp), intent(in)  dt,
real(rp), dimension(ka,ia,ja,vs:ve), intent(inout)  var 
)

Definition at line 270 of file scale_atmos_dyn_tinteg_rkcommon.F90.

270  implicit none
271 
272  type(RKInfo), intent(inout) :: this
273  integer, intent(in) :: io(this%var_num)
274  integer, intent(in) :: jo(this%var_num)
275  integer, intent(in) :: ko(this%var_num)
276  real(RP), intent(in) :: dt
277  integer, intent(in) :: vs, ve
278  real(RP), intent(inout) :: var(KA,IA,JA,vs:ve)
279 
280  integer :: i, j, k, iv, rks
281  real(RP) :: var0
282 
283  real(RP) :: b_(this%stage_num)
284  !--------------------------------------
285 
286  b_(:) = dt * this%rkcoef_b(:)
287 
288  !$omp parallel private(k,i,j,var0)
289 
290  rks = this%stage_num
291  do iv=vs, ve
292  !$omp do private(k,i,j,var0) collapse(2)
293  do j=js-jo(iv), je
294  do i=is-io(iv), ie
295  do k=ks-ko(iv), ke
296  var0 = this%work0(k,i,j,iv)
297  var(k,i,j,iv) = var0 + b_(rks) * ( this%work(k,i,j,iv,rks) - var0 )
298  end do
299  end do
300  end do
301  end do
302  do rks=1, this%stage_num-1
303  do iv=vs, ve
304  !$omp do private(k,i,j,var0) collapse(2)
305  do j=js-jo(iv), je
306  do i=is-io(iv), ie
307  do k=ks-ko(iv), ke
308  var(k,i,j,iv) = var(k,i,j,iv) + b_(rks) * this%work(k,i,j,iv,rks)
309  end do
310  end do
311  end do
312  end do
313  end do
314  !$omp end parallel
315 
316  return

References scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

◆ atmos_dyn_tinteg_rkcommon_updateflux()

subroutine, public scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_updateflux ( type(rkinfo), intent(inout)  this,
integer, intent(in)  nowstage,
integer, intent(in)  io,
integer, intent(in)  jo,
integer, intent(in)  ko,
integer, intent(in)  va_,
real(rp), dimension(ka,ia,ja,va_), intent(inout)  flux 
)

Definition at line 320 of file scale_atmos_dyn_tinteg_rkcommon.F90.

320  implicit none
321 
322  type(RKInfo), intent(inout) :: this
323  integer, intent(in) :: nowstage
324  integer, intent(in) :: va_
325  integer, intent(in) :: io, jo, ko
326  real(RP), intent(inout) :: flux(KA,IA,JA,va_)
327 
328  integer :: i, j, k, iv
329  !--------------------------------------
330 
331  if ( nowstage == 1) then
332  !$omp parallel do private(iv,k,j,i) collapse(3)
333  do iv=1, va_
334  do j=js-jo, je
335  do i=is-io, ie
336  do k=ks-ko, ke
337  flux(k,i,j,iv) = this%rkcoef_b(nowstage) * this%buf(k,i,j,iv)
338  end do
339  end do
340  end do
341  end do
342  else
343  !$omp parallel do private(iv,k,j,i) collapse(3)
344  do iv=1, va_
345  do j=js-jo, je
346  do i=is-io, ie
347  do k=ks-ko, ke
348  flux(k,i,j,iv) = flux(k,i,j,iv) + this%rkcoef_b(nowstage) * this%buf(k,i,j,iv)
349  end do
350  end do
351  end do
352  end do
353  end if
354 
355  return

References scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
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_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56