MONC
Data Types | Functions/Subroutines
stencil_mod Module Reference

Performs the interpolation between the primal and dual grids via a stencil approach. For performance reasons, for each field we store the entirety of the y dimension and the number of x slices required by the stencil. Therefore a new interpolation is simply the calculation of one point and reuse of existing computed points (unless this is the first x or y.) The applicable interpolation stenciled data is copied out, with :,1,1 being the central point, minus and plus in each y and x dimension as determined by the stencil size. This is done so that the stencil size can easily be different for each flow field and this is the case with u (where we need u-2 in the X.) More...

Data Types

type  grid_stencil_type
 Configuration for a specific stencil interpolation to perform. More...
 

Functions/Subroutines

type(grid_stencil_type) function, public create_stencil (local_grid, fields, nfields, interpolations_to_perform, sizes, xdim, ydim)
 Creates a stencil configuration which will then be used for interpolation. More...
 
subroutine, public free_stencil (stencil)
 Frees up the memory allocated to a stencil. More...
 
subroutine, public interpolate_to_dual (local_grid, field, stencil, x, y, interpolated_fields, interpolation_id)
 Interpolates the (vector) flow fields from the primal to dual grid based upon a specific field interpolation. More...
 
subroutine calculate_interpolated_cell_value (local_grid, j, i, interpolate_in_dimension, field, column_stencil_values)
 Calculates the interpolated values for each point in a column. The exact method of interpolation depends on whether we are interpolating in the w, y or x dimension. More...
 
subroutine copy_interpolated_to_prognostic (interpolated_fields, stencil, local_grid, x, y, interpolation_id)
 Copies stencil number of interpolated values from the store used by this functionality to preallocated prognostics. THe central point is always :,1,1 in the prognostic, with negative and in the x and y. More...
 

Detailed Description

Performs the interpolation between the primal and dual grids via a stencil approach. For performance reasons, for each field we store the entirety of the y dimension and the number of x slices required by the stencil. Therefore a new interpolation is simply the calculation of one point and reuse of existing computed points (unless this is the first x or y.) The applicable interpolation stenciled data is copied out, with :,1,1 being the central point, minus and plus in each y and x dimension as determined by the stencil size. This is done so that the stencil size can easily be different for each flow field and this is the case with u (where we need u-2 in the X.)

Function/Subroutine Documentation

◆ calculate_interpolated_cell_value()

subroutine stencil_mod::calculate_interpolated_cell_value ( type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  j,
integer, intent(in)  i,
logical, dimension(3)  interpolate_in_dimension,
type(prognostic_field_ptr_type), intent(inout)  field,
real(kind=default_precision), dimension(:), intent(inout)  column_stencil_values 
)
private

Calculates the interpolated values for each point in a column. The exact method of interpolation depends on whether we are interpolating in the w, y or x dimension.

Parameters
local_gridThe local grid
jCurrent local index in the y dimension
iCurrent local index in the x dimension
interpolate_in_dimensionWhich dimensions to interpolate in
fieldThe field to interpolate from
column_stencil_valuesThe result of the interpolation for the specific column

Definition at line 148 of file stencil.F90.

148  type(local_grid_type), intent(inout) :: local_grid
149  integer, intent(in) :: j,i
150  logical, dimension(3) :: interpolate_in_dimension
151  type(prognostic_field_ptr_type), intent(inout) :: field
152  real(kind=DEFAULT_PRECISION), dimension(:), intent(inout) :: column_stencil_values
153 
154  integer :: k, y_interpol, x_interpol
155 
156  y_interpol=merge(j+1, j, interpolate_in_dimension(2))
157  x_interpol=merge(i+1, i, interpolate_in_dimension(3))
158 
159  do k=1, local_grid%size(z_index)
160  if (y_interpol .lt. 1 .or. y_interpol .gt. local_grid%local_domain_end_index(y_index)+ local_grid%halo_size(y_index) &
161  .or. x_interpol .lt. 1 .or. x_interpol .gt. local_grid%local_domain_end_index(x_index)+ local_grid%halo_size(x_index) &
162  .or. (interpolate_in_dimension(1) .and. k == local_grid%size(z_index))) then
163  column_stencil_values(k)=field%ptr%data(k,j,i)
164  else
165  column_stencil_values(k)=0.5_default_precision*(&
166  field%ptr%data(k,j,i) + field%ptr%data(merge(k+1, k, interpolate_in_dimension(1)), y_interpol, x_interpol))
167  end if
168  end do
Here is the caller graph for this function:

◆ copy_interpolated_to_prognostic()

subroutine stencil_mod::copy_interpolated_to_prognostic ( type(prognostic_field_type), dimension(:), intent(inout), allocatable  interpolated_fields,
type(grid_stencil_type), intent(inout)  stencil,
type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  x,
integer, intent(in)  y,
integer, intent(in)  interpolation_id 
)
private

Copies stencil number of interpolated values from the store used by this functionality to preallocated prognostics. THe central point is always :,1,1 in the prognostic, with negative and in the x and y.

Parameters
interpolated_fieldsThe interpolated field prognostics that will be updated
stencilThe stencil configuration
local_gridLocal grid
xLocal x index of this column
yLocal y index of this column
interpolation_idThe ID of this interpolation which is used to retrieve information from the stencil configuration

Definition at line 180 of file stencil.F90.

180  type(prognostic_field_type), dimension(:), allocatable, intent(inout) :: interpolated_fields
181  type(local_grid_type), intent(inout) :: local_grid
182  type(grid_stencil_type), intent(inout) :: stencil
183  integer, intent(in) :: x, y, interpolation_id
184 
185  integer :: i, y_min_src, y_max_src, x_max_src, y_min_tgt, y_max_tgt, x_min_tgt, x_max_tgt
186 
187  do i=1,stencil%nfields
188  y_min_src=y-stencil%sizes(i,2)
189  y_min_tgt=1-stencil%sizes(i,2)
190  y_max_src=y+stencil%sizes(i,2)
191  y_max_tgt=1+stencil%sizes(i,2)
192  x_min_tgt=1-stencil%sizes(i,1)
193  x_max_src=stencil%sizes(i,1)*2+1
194  x_max_tgt=1+stencil%sizes(i,1)
195  if (y_min_src .lt. 1) then
196  y_min_tgt=y_min_tgt+(1-y_min_src)
197  y_min_src=1
198  end if
199 
200  interpolated_fields(i)%data(:,y_min_tgt:y_max_tgt, x_min_tgt:x_max_tgt)=&
201  stencil%interpolated_fields(interpolation_id, i)%data(:, y_min_src:y_max_src, 1:x_max_src)
202  end do
Here is the caller graph for this function:

◆ create_stencil()

type(grid_stencil_type) function, public stencil_mod::create_stencil ( type(local_grid_type), intent(inout)  local_grid,
type(prognostic_field_ptr_type), dimension(:), intent(inout)  fields,
integer, intent(in)  nfields,
integer, intent(in)  interpolations_to_perform,
integer, dimension(:,:), intent(in)  sizes,
logical, intent(in)  xdim,
logical, intent(in)  ydim 
)

Creates a stencil configuration which will then be used for interpolation.

Parameters
local_gridThe local grid
fieldsThe fields to perform the interpolation upon
nfieldsThe number of fields included in a specific interpolation
interpolations_to_performThe number of distinct interpolations that will be performed
sizesThe size of the stencil, for each field, in X and Y dimensions
xdimWhether the X dimension is active
ydimWhether the Y dimension is active

Definition at line 37 of file stencil.F90.

37  type(local_grid_type), intent(inout) :: local_grid
38  type(prognostic_field_ptr_type), dimension(:), intent(inout) :: fields
39  integer, dimension(:,:), intent(in) :: sizes
40  integer, intent(in) :: nfields, interpolations_to_perform
41  logical, intent(in) :: xdim, ydim
42 
43  integer :: j,i
44 
45  allocate(create_stencil%fields(0:nfields), create_stencil%interpolated_fields(interpolations_to_perform,nfields))
46  allocate(create_stencil%sizes(nfields,2))
47  do j=1,nfields
48  create_stencil%fields(j) = fields(j)
49  create_stencil%sizes(j,1) = merge(sizes(j,1), 0, xdim)
50  create_stencil%sizes(j,2) = merge(sizes(j,2), 0, ydim)
51  do i=1,interpolations_to_perform
52  allocate(create_stencil%interpolated_fields(i, j)%data(local_grid%size(z_index), &
53  local_grid%size(y_index)+4, create_stencil%sizes(j,1)*2+1))
54  end do
55  end do
56  create_stencil%nfields=nfields
57  create_stencil%max_y_point=(local_grid%size(y_index)+local_grid%halo_size(y_index) * 2) - 1
58  create_stencil%max_x_point=(local_grid%size(x_index)+local_grid%halo_size(x_index) * 2) - 1
Here is the caller graph for this function:

◆ free_stencil()

subroutine, public stencil_mod::free_stencil ( type(grid_stencil_type), intent(inout)  stencil)

Frees up the memory allocated to a stencil.

Parameters
stencilThe stencil that we will free up

Definition at line 64 of file stencil.F90.

64  type(grid_stencil_type), intent(inout) :: stencil
65 
66  integer :: j, i
67 
68  if (allocated(stencil%interpolated_fields)) then
69  do i=1, size(stencil%interpolated_fields,1)
70  do j=1, size(stencil%interpolated_fields,2)
71  if (allocated(stencil%interpolated_fields(i,j)%data)) deallocate(stencil%interpolated_fields(i,j)%data)
72  end do
73  end do
74  deallocate(stencil%fields, stencil%interpolated_fields)
75  end if
76  if (allocated(stencil%sizes)) deallocate(stencil%sizes)

◆ interpolate_to_dual()

subroutine, public stencil_mod::interpolate_to_dual ( type(local_grid_type), intent(inout)  local_grid,
type(prognostic_field_type), intent(inout)  field,
type(grid_stencil_type), intent(inout)  stencil,
integer, intent(in)  x,
integer, intent(in)  y,
type(prognostic_field_type), dimension(:), intent(inout), allocatable  interpolated_fields,
integer, intent(in)  interpolation_id 
)

Interpolates the (vector) flow fields from the primal to dual grid based upon a specific field interpolation.

Parameters
local_gridThe local grid
fieldThe field to base this interpolation upon
stencilStencil configuration
xLocal column index in X
yLocal column index in Y
interpolated_fieldsInterpolated prognostics that this will write the output into, central (x,y) point is 1,1
interpolation_idId of the interpolation

Definition at line 88 of file stencil.F90.

88  type(prognostic_field_type), intent(inout) :: field
89  type(local_grid_type), intent(inout) :: local_grid
90  type(grid_stencil_type), intent(inout) :: stencil
91  type(prognostic_field_type), dimension(:), allocatable, intent(inout) :: interpolated_fields
92  integer, intent(in) :: x, y, interpolation_id
93 
94  integer :: i,j,n, x_top, base_x, y_start, y_top, x_start, abs_x_top
95  logical, dimension(3) :: interpolate_in_dimension
96 
97  interpolate_in_dimension = get_field_interpolation_index(field)
98 
99  do n=1,stencil%nfields
100  base_x=x-(stencil%sizes(n,1)+1)
101  abs_x_top=x+stencil%sizes(n,1)
102  x_start=1
103  if (base_x .lt. 0) x_start=x_start+(0-base_x)
104  x_top=stencil%sizes(n,1)*2+1
105  y_start=y-stencil%sizes(n,2)
106  y_top=y+stencil%sizes(n,2)
107  if (y_start .lt. 1) y_start=1
108 
109  if (x==2) then
110  ! If this is the first x slice then compute all points on y dimension up to x top - 1
111  do i=x_start,x_top-1
112  do j=y_start, y_top
113  call calculate_interpolated_cell_value(local_grid, j, base_x+i, interpolate_in_dimension, stencil%fields(n), &
114  stencil%interpolated_fields(interpolation_id,n)%data(:,j,i))
115  end do
116  end do
117  end if
118 
119  if (y==2) then
120  ! If this is the first column in y then shift X down (as only store required data) and calculate y to y top -1
121  if (x/=2) then
122  do i=1, x_top-1
123  stencil%interpolated_fields(interpolation_id,n)%data(:,:,i)=&
124  stencil%interpolated_fields(interpolation_id,n)%data(:,:,i+1)
125  end do
126  end if
127  do j=y_start, y_top-1
128  call calculate_interpolated_cell_value(local_grid, j, abs_x_top, interpolate_in_dimension, stencil%fields(n), &
129  stencil%interpolated_fields(interpolation_id,n)%data(:,j,x_top))
130  end do
131  end if
132  ! Calculate y top, x top which is done for each column
133  call calculate_interpolated_cell_value(local_grid, y_top, abs_x_top, interpolate_in_dimension, stencil%fields(n), &
134  stencil%interpolated_fields(interpolation_id,n)%data(:,y_top,x_top))
135  end do
136  call copy_interpolated_to_prognostic(interpolated_fields, stencil, local_grid, x, y, interpolation_id)
Here is the call graph for this function: