MONC
lateral_bcs.F90
Go to the documentation of this file.
1 
6  use state_mod, only : model_state_type
7  use grids_mod, only : x_index, y_index, z_index
10  implicit none
11 
12 #ifndef TEST_MODE
13  private
14 #endif
15 
16  INTEGER :: lbc_type_x, lbc_type_y ! probably need to add these to state
17 
18  INTEGER :: idepth, jdepth ! depth of lateral bcs in i and j directions
19 
20  INTEGER, PARAMETER :: lbc_periodic = 1
21  INTEGER, PARAMETER :: lbc_rigid = 2
22  INTEGER, PARAMETER :: lbc_blend = 3
23 
25 
26 contains
27 
31  lateral_bcs_get_descriptor%name="lateral_bcs"
32  lateral_bcs_get_descriptor%version=0.1
35  end function lateral_bcs_get_descriptor
36 
39  subroutine initialisation_callback(current_state)
40  type(model_state_type), target, intent(inout) :: current_state
41 
42  ! Hard wiring these for now
43  idepth=3
44  jdepth=2
45 
48 
49  end subroutine initialisation_callback
50 
54  subroutine timestep_callback(current_state)
55  type(model_state_type), target, intent(inout) :: current_state
56 
57  INTEGER :: istart, iend, jstart, jend, kstart, kend
58 
59  INTEGER :: i,j,k,iq
60 
61  LOGICAL :: is_east, is_west, is_north, is_south
62 
63  if (current_state%timestep==1) return
64 
65  is_west=current_state%parallel%wrapped_around(x_index,1)
66  is_east=current_state%parallel%wrapped_around(x_index,2)
67  is_south=current_state%parallel%wrapped_around(y_index,1)
68  is_north=current_state%parallel%wrapped_around(y_index,2)
69 
70  ! k points
71  if (is_west .or. is_east .or. is_north .or. is_south)then
72  kstart=current_state%local_grid%local_domain_start_index(z_index)
73  kend=current_state%local_grid%local_domain_end_index(z_index)
74  end if
75 
76  ! Eastern boundary
77  if (is_east)then
78 
79  iend=current_state%local_grid%local_domain_end_index(x_index)
80  istart = iend - idepth + 1
81  jstart=current_state%local_grid%local_domain_start_index(y_index)
82  jend=current_state%local_grid%local_domain_end_index(y_index)
83 
84  if (lbc_type_x==lbc_periodic)then
85  ! how should this be done best Luis?
86  else if (lbc_type_x==lbc_rigid)then
87  ! No communications needed assuming idepth and jdepth fit
88  ! on a single pe (what if they don't?)
89  call apply_rigid(current_state,istart,iend,jstart,jend,kstart,kend)
90  ! propagate winds to halos
91  do i=iend+1,iend+current_state%local_grid%halo_size(x_index)
92  current_state%u%data(:,:,i) = current_state%u%data(:,:,iend)
93  current_state%zu%data(:,:,i) = current_state%zu%data(:,:,iend)
94  current_state%w%data(:,:,i) = current_state%w%data(:,:,iend)
95  current_state%zw%data(:,:,i) = current_state%zw%data(:,:,iend)
96  end do
97 
98  else if (lbc_type_x==lbc_blend)then
99  ! To be coded....
100  end if
101 
102  end if
103 
104  ! Western boundary
105  if (is_west)then
106 
107  istart=current_state%local_grid%local_domain_start_index(x_index)
108  iend=istart + idepth - 1
109  jstart=current_state%local_grid%local_domain_start_index(y_index)
110  jend=current_state%local_grid%local_domain_end_index(y_index)
111 
112 
113  if (lbc_type_x==lbc_periodic)then
114  ! how should this be done best Luis?
115  else if (lbc_type_x==lbc_rigid)then
116  ! No communications needed assuming idepth and jdepth fit
117  ! on a single pe (what if they don't?)
118  call apply_rigid(current_state,istart,iend,jstart,jend,kstart,kend)
119  ! propagate winds to halos
120  do i=istart-current_state%local_grid%halo_size(x_index),istart-1
121  current_state%u%data(:,:,i) = current_state%u%data(:,:,istart)
122  current_state%zu%data(:,:,i) = current_state%zu%data(:,:,istart)
123  current_state%w%data(:,:,i) = current_state%w%data(:,:,istart)
124  current_state%zw%data(:,:,i) = current_state%zw%data(:,:,istart)
125  end do
126  else if (lbc_type_x==lbc_blend)then
127  ! To be coded....
128  end if
129  end if
130 
131 
132  ! Northern boundary
133  if (is_north)then
134 
135  istart=current_state%local_grid%local_domain_start_index(x_index)
136  iend=current_state%local_grid%local_domain_end_index(x_index)
137  jend=current_state%local_grid%local_domain_end_index(y_index)
138  jstart=jend - jdepth + 1
139 
140  if (lbc_type_y==lbc_periodic)then
141  ! how should this be done best Luis?
142  else if (lbc_type_y==lbc_rigid)then
143  ! No communications needed assuming idepth and jdepth fit
144  ! on a single pe (what if they don't?)
145  call apply_rigid(current_state,istart,iend,jstart,jend,kstart,kend)
146  do j=jend+1,jend+current_state%local_grid%halo_size(y_index)
147  current_state%v%data(:,:,j) = current_state%v%data(:,:,jend)
148  current_state%zv%data(:,:,j) = current_state%zv%data(:,:,jend)
149  current_state%w%data(:,:,j) = current_state%w%data(:,:,jend)
150  current_state%zw%data(:,:,j) = current_state%zw%data(:,:,jend)
151  end do
152  else if (lbc_type_y==lbc_blend)then
153  ! To be coded....
154  end if
155 
156  end if
157 
158 
159  ! Southern boundary
160  if (is_south)then
161 
162  istart=current_state%local_grid%local_domain_start_index(x_index)
163  iend=current_state%local_grid%local_domain_end_index(x_index)
164  jstart=current_state%local_grid%local_domain_start_index(y_index)
165  jend=jstart + jdepth - 1
166 
167  if (lbc_type_y==lbc_periodic)then
168  ! how should this be done best Luis?
169  else if (lbc_type_y==lbc_rigid)then
170  ! No communications needed assuming idepth and jdepth fit
171  ! on a single pe (what if they don't?)
172  call apply_rigid(current_state,istart,iend,jstart,jend,kstart,kend)
173  do j=jstart-current_state%local_grid%halo_size(y_index),jstart-1
174  current_state%v%data(:,:,j) = current_state%v%data(:,:,jstart)
175  current_state%zv%data(:,:,j) = current_state%zv%data(:,:,jstart)
176  current_state%w%data(:,:,j) = current_state%w%data(:,:,jstart)
177  current_state%zw%data(:,:,j) = current_state%zw%data(:,:,jstart)
178  end do
179  else if (lbc_type_y==lbc_blend)then
180  ! To be coded....
181  end if
182 
183  end if
184  end subroutine timestep_callback
185 
186  subroutine apply_rigid(current_state,istart,iend,jstart,jend,kstart,kend)
187  type(model_state_type), target, intent(inout) :: current_state
188 
189  INTEGER :: istart, iend, jstart, jend, kstart, kend
190 
191  INTEGER :: i,j,k,iq
192 
193  do i=istart, iend
194  do j=jstart, jend
195  do k=kstart,kend
196  ! Update to old values, i.e. fixed in time
197 #ifdef U_ACTIVE
198 ! current_state%u%data(k,j,i) = current_state%zu%data(k,j,i)
199 #endif
200 #ifdef V_ACTIVE
201 ! current_state%v%data(k,j,i) = current_state%zv%data(k,j,i)
202 #endif
203 #ifdef W_ACTIVE
204 ! current_state%w%data(k,j,i) = current_state%zw%data(k,j,i)
205 #endif
206  current_state%th%data(k,j,i) = current_state%zth%data(k,j,i)
207  end do
208  end do
209  end do
210  do iq = 1, current_state%number_q_fields
211  do i=istart, iend
212  do j=jstart, jend
213  do k=kstart,kend
214  current_state%q(iq)%data(k,j,i) = current_state%zq(iq)%data(k,j,i)
215  end do
216  end do
217  end do
218  end do
219 
220  end subroutine apply_rigid
221 
222 end module lateral_bcs_mod
lateral_bcs_mod
Applies lateral boundary conditions. Note that these boundary conditions only make sense if bi-period...
Definition: lateral_bcs.F90:4
lateral_bcs_mod::timestep_callback
subroutine timestep_callback(current_state)
Timestep callback which applies lateral boundary conditions.
Definition: lateral_bcs.F90:55
lateral_bcs_mod::lbc_periodic
integer, parameter lbc_periodic
Definition: lateral_bcs.F90:20
lateral_bcs_mod::initialisation_callback
subroutine initialisation_callback(current_state)
Initialisation callback. Sets up boundary depth and other configurable options.
Definition: lateral_bcs.F90:40
lateral_bcs_mod::lbc_type_x
integer lbc_type_x
Definition: lateral_bcs.F90:16
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.F90:14
lateral_bcs_mod::lbc_blend
integer, parameter lbc_blend
Definition: lateral_bcs.F90:22
monc_component_mod
Interfaces and types that MONC components must specify.
Definition: monc_component.F90:6
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
state_mod::model_state_type
The ModelState which represents the current state of a run.
Definition: state.F90:39
lateral_bcs_mod::lateral_bcs_get_descriptor
type(component_descriptor_type) function, public lateral_bcs_get_descriptor()
Provides a description of this component for the core to register.
Definition: lateral_bcs.F90:31
lateral_bcs_mod::idepth
integer idepth
Definition: lateral_bcs.F90:18
logging_mod
Logging utility.
Definition: logging.F90:2
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
lateral_bcs_mod::apply_rigid
subroutine apply_rigid(current_state, istart, iend, jstart, jend, kstart, kend)
Definition: lateral_bcs.F90:187
lateral_bcs_mod::jdepth
integer jdepth
Definition: lateral_bcs.F90:18
lateral_bcs_mod::lbc_type_y
integer lbc_type_y
Definition: lateral_bcs.F90:16
grids_mod
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
monc_component_mod::component_descriptor_type
Description of a component.
Definition: monc_component.F90:42
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
state_mod
The model state which represents the current state of a run.
Definition: state.F90:2
lateral_bcs_mod::lbc_rigid
integer, parameter lbc_rigid
Definition: lateral_bcs.F90:21