FMS  2025.04
Flexible Modeling System
coupler_types.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup coupler_types_mod coupler_types_mod
19 !> @ingroup coupler
20 !> @brief This module contains type declarations for the coupler.
21 !> @author Richard Slater, John Dunne
22 
23 !> @addtogroup coupler_types_mod
24 !> @{
25 
26 !> @addtogroup coupler_types_mod
27 !> @{
28 
29  !> @brief Rescales the fields in the fields in the elements of a coupler_2d_bc_type
30  !!
31  !! Rescales the fields in the elements of a coupler_2d_bc_type by multiplying by a factor scale.
32  !! If scale is 0, this is a direct assignment to 0, so that NaNs will not persist.
33  subroutine ct_rescale_data_2d_(var, scale, halo_size, bc_index, field_index,&
34  & exclude_flux_type, only_flux_type, pass_through_ice)
35  type(coupler_2d_bc_type), intent(inout) :: var !< The BC_type structure whose fields are being rescaled
36  real(FMS_CP_KIND_), intent(in) :: scale !< A scaling factor to multiply fields by
37  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default or
38  !! the full arrays if scale is 0.
39  integer, optional, intent(in) :: bc_index !< The index of the boundary condition
40  !! that is being copied
41  integer, optional, intent(in) :: field_index !< The index of the field in the
42  !! boundary condition that is being copied
43  character(len=*), optional, intent(in) :: exclude_flux_type !< A string describing which types
44  !! of fluxes to exclude from this copy.
45  character(len=*), optional, intent(in) :: only_flux_type !< A string describing which types
46  !! of fluxes to include from this copy.
47  logical, optional, intent(in) :: pass_through_ice !< If true, only copy BCs whose
48  !! value of pass_through ice matches this
49 
50  logical :: do_bc
51  integer :: i, j, m, n, n1, n2, halo
52  integer, parameter :: kindl = fms_cp_kind_
53 
54  if (present(bc_index)) then
55  if (bc_index > var%num_bcs)&
56  & call mpp_error(fatal, "CT_rescale_data_2d: bc_index is present and exceeds var%num_bcs.")
57  if (present(field_index)) then ; if (field_index > var%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
58  & call mpp_error(fatal, "CT_rescale_data_2d: field_index is present and exceeds num_fields for" //&
59  & trim(var%FMS_CP_BC_TYPE_(bc_index)%name) )
60  endif
61  elseif (present(field_index)) then
62  call mpp_error(fatal, "CT_rescale_data_2d: bc_index must be present if field_index is present.")
63  endif
64 
65  halo = 0
66  if (present(halo_size)) halo = halo_size
67 
68  n1 = 1
69  n2 = var%num_bcs
70  if (present(bc_index)) then
71  n1 = bc_index
72  n2 = bc_index
73  endif
74 
75  if (n2 >= n1) then
76  ! A more consciencious implementation would include a more descriptive error messages.
77  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
78  & call mpp_error(fatal, "CT_rescale_data_2d: Excessive i-direction halo size.")
79  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
80  & call mpp_error(fatal, "CT_rescale_data_2d: Excessive j-direction halo size.")
81  endif
82 
83  do n = n1, n2
84  do_bc = .true.
85  if (do_bc .and. present(exclude_flux_type))&
86  & do_bc = .not.(trim(var%FMS_CP_BC_TYPE_(n)%flux_type) == trim(exclude_flux_type))
87  if (do_bc .and. present(only_flux_type))&
88  & do_bc = (trim(var%FMS_CP_BC_TYPE_(n)%flux_type) == trim(only_flux_type))
89  if (do_bc .and. present(pass_through_ice))&
90  & do_bc = (pass_through_ice .eqv. var%FMS_CP_BC_TYPE_(n)%pass_through_ice)
91  if (.not.do_bc) cycle
92 
93  do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
94  if (present(field_index)) then
95  if (m /= field_index) cycle
96  endif
97  if ( associated(var%FMS_CP_BC_TYPE_(n)%field(m)%values) ) then
98  if (scale == 0.0_kindl) then
99  if (present(halo_size)) then
100  do j=var%jsc-halo,var%jec+halo
101  do i=var%isc-halo,var%iec+halo
102  var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j) = 0.0_kindl
103  enddo
104  enddo
105  else
106  var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:) = 0.0_kindl
107  endif
108  else
109  do j=var%jsc-halo,var%jec+halo
110  do i=var%isc-halo,var%iec+halo
111  var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j) = scale * var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j)
112  enddo
113  enddo
114  endif
115  endif
116  enddo
117  enddo
118  end subroutine ct_rescale_data_2d_
119 
120  !> @brief Rescales the fields in the elements of a coupler_3d_bc_type
121  !!
122  !! This subroutine rescales the fields in the elements of a coupler_3d_bc_type by multiplying by a
123  !! factor scale. If scale is 0, this is a direct assignment to 0, so that NaNs will not persist.
124  subroutine ct_rescale_data_3d_(var, scale, halo_size, bc_index, field_index,&
125  & exclude_flux_type, only_flux_type, pass_through_ice)
126  type(coupler_3d_bc_type), intent(inout) :: var !< The BC_type structure whose fields are being rescaled
127  real(FMS_CP_KIND_), intent(in) :: scale !< A scaling factor to multiply fields by
128  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default or
129  !! the full arrays if scale is 0.
130  integer, optional, intent(in) :: bc_index !< The index of the boundary condition
131  !! that is being copied
132  integer, optional, intent(in) :: field_index !< The index of the field in the
133  !! boundary condition that is being copied
134  character(len=*), optional, intent(in) :: exclude_flux_type !< A string describing which types
135  !! of fluxes to exclude from this copy.
136  character(len=*), optional, intent(in) :: only_flux_type !< A string describing which types of
137  !! fluxes to include from this copy.
138  logical, optional, intent(in) :: pass_through_ice !< If true, only copy BCs whose
139  !! value of pass_through ice matches this
140 
141  logical :: do_bc
142  integer :: i, j, k, m, n, n1, n2, halo
143  integer, parameter :: kindl = fms_cp_kind_
144 
145  if (present(bc_index)) then
146  if (bc_index > var%num_bcs)&
147  & call mpp_error(fatal, "CT_rescale_data_2d: bc_index is present and exceeds var%num_bcs.")
148  if (present(field_index)) then ; if (field_index > var%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
149  & call mpp_error(fatal, "CT_rescale_data_2d: field_index is present and exceeds num_fields for" //&
150  & trim(var%FMS_CP_BC_TYPE_(bc_index)%name) )
151  endif
152  elseif (present(field_index)) then
153  call mpp_error(fatal, "CT_rescale_data_2d: bc_index must be present if field_index is present.")
154  endif
155 
156  halo = 0
157  if (present(halo_size)) halo = halo_size
158 
159  n1 = 1
160  n2 = var%num_bcs
161  if (present(bc_index)) then
162  n1 = bc_index
163  n2 = bc_index
164  endif
165 
166  if (n2 >= n1) then
167  ! A more consciencious implementation would include a more descriptive error messages.
168  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
169  & call mpp_error(fatal, "CT_rescale_data_3d: Excessive i-direction halo size.")
170  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
171  & call mpp_error(fatal, "CT_rescale_data_3d: Excessive j-direction halo size.")
172  endif
173 
174  do n = n1, n2
175  do_bc = .true.
176  if (do_bc .and. present(exclude_flux_type))&
177  & do_bc = .not.(trim(var%FMS_CP_BC_TYPE_(n)%flux_type) == trim(exclude_flux_type))
178  if (do_bc .and. present(only_flux_type))&
179  & do_bc = (trim(var%FMS_CP_BC_TYPE_(n)%flux_type) == trim(only_flux_type))
180  if (do_bc .and. present(pass_through_ice))&
181  & do_bc = (pass_through_ice .eqv. var%FMS_CP_BC_TYPE_(n)%pass_through_ice)
182  if (.not.do_bc) cycle
183 
184  do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
185  if (present(field_index)) then
186  if (m /= field_index) cycle
187  endif
188  if ( associated(var%FMS_CP_BC_TYPE_(n)%field(m)%values) ) then
189  if (scale == 0.0_kindl) then
190  if (present(halo_size)) then
191  do k=var%ks,var%ke
192  do j=var%jsc-halo,var%jec+halo
193  do i=var%isc-halo,var%iec+halo
194  var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j,k) = 0.0_kindl
195  enddo
196  enddo
197  enddo
198  else
199  var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:,:) = 0.0_kindl
200  endif
201  else
202  do k=var%ks,var%ke
203  do j=var%jsc-halo,var%jec+halo
204  do i=var%isc-halo,var%iec+halo
205  var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j,k) = scale * var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j,k)
206  enddo
207  enddo
208  enddo
209  endif
210  endif
211  enddo
212  enddo
213  end subroutine ct_rescale_data_3d_
214 
215  !> @brief Increment data in the elements of a coupler_2d_bc_type with weighted averages of elements of a
216  !! coupler_3d_bc_type
217  !!
218  !! Increments the data in the elements of a coupler_2d_bc_type with the weighed average of the
219  !! elements of a coupler_3d_bc_type. Both must have the same horizontal array sizes and the
220  !! normalized weight array must match the array sizes of the coupler_3d_bc_type.
221  !!
222  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
223  !! @throw FATAL, "field_index is present and exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
224  !! @throw FATAL, "bc_index must be present if field_index is present."
225  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
226  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
227  !! @throw FATAL, "There is an k-direction computational domain size mismatch."
228  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
229  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
230  !! @throw FATAL, "weights array must be the i-size of a computational or data domain."
231  !! @throw FATAL, "weights array must be the j-size of a computational or data domain."
232  subroutine ct_increment_data_2d_3d_(var_in, weights, var, halo_size, bc_index, field_index,&
233  & scale_factor, scale_prev, exclude_flux_type, only_flux_type, pass_through_ice)
234  type(coupler_3d_bc_type), intent(in) :: var_in !< BC_type structure with the data to add to the other type
235  real(FMS_CP_KIND_), dimension(:,:,:), intent(in) :: weights !< An array of normalized weights for the 3d-data
236  !! to increment the 2d-data. There is no renormalization,
237  !! so if the weights do not sum to 1 in the 3rd dimension
238  !! there may be adverse consequences!
239  type(coupler_2d_bc_type), intent(inout) :: var !< The BC_type structure whose fields are being incremented
240  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
241  integer, optional, intent(in) :: bc_index !< The index of the boundary condition
242  !! that is being copied
243  integer, optional, intent(in) :: field_index !< The index of the field in the
244  !! boundary condition that is being copied
245  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that
246  !! is being added
247  real(FMS_CP_KIND_), optional, intent(in) :: scale_prev !< A scaling factor for the data that
248  !! is already here and is being added
249  character(len=*), optional, intent(in) :: exclude_flux_type !< A string describing which types
250  !! of fluxes to exclude from this increment.
251  character(len=*), optional, intent(in) :: only_flux_type !< A string describing which types
252  !! of fluxes to include from this increment.
253  logical, optional, intent(in) :: pass_through_ice !< If true, only increment BCs whose
254  !! value of pass_through ice matches this
255 
256  real(FMS_CP_KIND_) :: scale, sc_prev
257  logical :: increment_bc
258  integer :: i, j, k, m, n, n1, n2, halo
259  integer :: io1, jo1, iow, jow, kow ! Offsets to account for different index conventions.
260  integer, parameter :: kindl = fms_cp_kind_
261 
262  scale = 1.0_kindl
263  if (present(scale_factor)) scale = scale_factor
264  sc_prev = 1.0_kindl
265  if (present(scale_prev)) sc_prev = scale_prev
266 
267  if (present(bc_index)) then
268  if (bc_index > var_in%num_bcs)&
269  & call mpp_error(fatal, "CT_increment_data_2d_3d: bc_index is present and exceeds var_in%num_bcs.")
270  if (present(field_index)) then ; if (field_index > var_in%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
271  & call mpp_error(fatal, "CT_increment_data_2d_3d: field_index is present and exceeds num_fields for" //&
272  & trim(var_in%FMS_CP_BC_TYPE_(bc_index)%name) )
273  endif
274  elseif (present(field_index)) then
275  call mpp_error(fatal, "CT_increment_data_2d_3d: bc_index must be present if field_index is present.")
276  endif
277 
278  halo = 0
279  if (present(halo_size)) halo = halo_size
280 
281  n1 = 1
282  n2 = var_in%num_bcs
283  if (present(bc_index)) then
284  n1 = bc_index
285  n2 = bc_index
286  endif
287 
288  if (n2 >= n1) then
289  ! A more consciencious implementation would include a more descriptive error messages.
290  if ((var_in%iec-var_in%isc) /= (var%iec-var%isc))&
291  & call mpp_error(fatal, &
292  & "CT_increment_data_2d_3d: There is an i-direction computational domain size mismatch.")
293  if ((var_in%jec-var_in%jsc) /= (var%jec-var%jsc))&
294  & call mpp_error(fatal, &
295  & "CT_increment_data_2d_3d: There is a j-direction computational domain size mismatch.")
296  if ((1+var_in%ke-var_in%ks) /= size(weights,3))&
297  & call mpp_error(fatal, &
298  & "CT_increment_data_2d_3d: There is a k-direction size mismatch with the weights array.")
299  if ((var_in%isc-var_in%isd < halo) .or. (var_in%ied-var_in%iec < halo))&
300  & call mpp_error(fatal, "CT_increment_data_2d_3d: Excessive i-direction halo size for the input structure.")
301  if ((var_in%jsc-var_in%jsd < halo) .or. (var_in%jed-var_in%jec < halo))&
302  & call mpp_error(fatal, "CT_increment_data_2d_3d: Excessive j-direction halo size for the input structure.")
303  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
304  & call mpp_error(fatal, "CT_increment_data_2d_3d: Excessive i-direction halo size for the output structure.")
305  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
306  & call mpp_error(fatal, "CT_increment_data_2d_3d: Excessive j-direction halo size for the output structure.")
307 
308  if ((1+var%iec-var%isc) == size(weights,1)) then
309  iow = 1 - var%isc
310  elseif ((1+var%ied-var%isd) == size(weights,1)) then
311  iow = 1 - var%isd
312  elseif ((1+var_in%ied-var_in%isd) == size(weights,1)) then
313  iow = 1 + (var_in%isc - var_in%isd) - var%isc
314  else
315  call mpp_error(fatal, &
316  & "CT_increment_data_2d_3d: weights array must be the i-size of a computational or data domain.")
317  endif
318  if ((1+var%jec-var%jsc) == size(weights,2)) then
319  jow = 1 - var%jsc
320  elseif ((1+var%jed-var%jsd) == size(weights,2)) then
321  jow = 1 - var%jsd
322  elseif ((1+var_in%jed-var_in%jsd) == size(weights,2)) then
323  jow = 1 + (var_in%jsc - var_in%jsd) - var%jsc
324  else
325  call mpp_error(fatal, &
326  & "CT_increment_data_2d_3d: weights array must be the j-size of a computational or data domain.")
327  endif
328 
329  io1 = var_in%isc - var%isc
330  jo1 = var_in%jsc - var%jsc
331  kow = 1 - var_in%ks
332  endif
333 
334  do n = n1, n2
335  increment_bc = .true.
336  if (increment_bc .and. present(exclude_flux_type))&
337  & increment_bc = .not.(trim(var_in%FMS_CP_BC_TYPE_(n)%flux_type) == trim(exclude_flux_type))
338  if (increment_bc .and. present(only_flux_type))&
339  & increment_bc = (trim(var_in%FMS_CP_BC_TYPE_(n)%flux_type) == trim(only_flux_type))
340  if (increment_bc .and. present(pass_through_ice))&
341  & increment_bc = (pass_through_ice .eqv. var_in%FMS_CP_BC_TYPE_(n)%pass_through_ice)
342  if (.not.increment_bc) cycle
343 
344  do m = 1, var_in%FMS_CP_BC_TYPE_(n)%num_fields
345  if (present(field_index)) then
346  if (m /= field_index) cycle
347  endif
348  if ( associated(var%FMS_CP_BC_TYPE_(n)%field(m)%values) ) then
349  do k=var_in%ks,var_in%ke
350  do j=var%jsc-halo,var%jec+halo
351  do i=var%isc-halo,var%iec+halo
352  var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j) = sc_prev * var%FMS_CP_BC_TYPE_(n)%field(m)%values(i,j) +&
353  & (scale * weights(i+iow,j+jow,k+kow)) * var_in%FMS_CP_BC_TYPE_(n)%field(m)%values(i+io1,j+io1,k)
354  enddo
355  enddo
356  enddo
357  endif
358  enddo
359  enddo
360  end subroutine ct_increment_data_2d_3d_
361 
362  !> @brief Extract a 2d field from a coupler_2d_bc_type
363  !!
364  !! Extract a single 2-d field from a coupler_2d_bc_type into a two-dimensional array.
365  !!
366  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
367  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
368  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
369  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
370  !! @throw FATAL, "Disordered i-dimension index bound list"
371  !! @throw FATAL, "Disordered j-dimension index bound list"
372  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
373  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
374  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
375  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
376  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
377  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
378  subroutine ct_extract_data_2d_(var_in, bc_index, field_index, array_out,&
379  & scale_factor, halo_size, idim, jdim)
380  type(coupler_2d_bc_type), intent(in) :: var_in !< BC_type structure with the data to extract
381  integer, intent(in) :: bc_index !< The index of the boundary condition
382  !! that is being copied
383  integer, intent(in) :: field_index !< The index of the field in the
384  !! boundary condition that is being copied
385  real(FMS_CP_KIND_), dimension(1:,1:), intent(out) :: array_out !< The recipient array for the field; its size
386  !! must match the size of the data being copied
387  !! unless idim and jdim are supplied.
388  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
389  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
390  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
391  !! the first dimension of the output array
392  !! in a non-decreasing list
393  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
394  !! the second dimension of the output array
395  !! in a non-decreasing list
396 
397  character(len=*), parameter :: error_header =&
398  & '==>Error from coupler_types_mod (CT_extract_data_2d):'
399  character(len=400) :: error_msg
400 
401  real(FMS_CP_KIND_) :: scale
402  integer :: i, j, halo, i_off, j_off
403  integer, parameter :: kindl = fms_cp_kind_
404 
405  if (bc_index <= 0) then
406  array_out(:,:) = 0.0_kindl
407  return
408  endif
409 
410  halo = 0
411  if (present(halo_size)) halo = halo_size
412  scale = 1.0_kindl
413  if (present(scale_factor)) scale = scale_factor
414 
415  if ((var_in%isc-var_in%isd < halo) .or. (var_in%ied-var_in%iec < halo))&
416  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
417  if ((var_in%jsc-var_in%jsd < halo) .or. (var_in%jed-var_in%jec < halo))&
418  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
419 
420  if (bc_index > var_in%num_bcs)&
421  & call mpp_error(fatal, trim(error_header)//" bc_index exceeds var_in%num_bcs.")
422  if (field_index > var_in%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
423  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
424  & trim(var_in%FMS_CP_BC_TYPE_(bc_index)%name) )
425 
426  ! Do error checking on the i-dimension and determine the array offsets.
427  if (present(idim)) then
428  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
429  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
430  call mpp_error(fatal, trim(error_msg))
431  endif
432  if (size(array_out,1) /= (1+idim(4)-idim(1))) then
433  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
434  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_out,1)
435  call mpp_error(fatal, trim(error_msg))
436  endif
437  if ((var_in%iec-var_in%isc) /= (idim(3)-idim(2)))&
438  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
439  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
440  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
441  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
442  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
443  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
444  & (2*halo + 1 + var_in%iec - var_in%isc)
445  call mpp_error(fatal, trim(error_msg))
446  endif
447 
448  i_off = (1-idim(1)) + (idim(2)-var_in%isc)
449  else
450  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
451  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
452  & size(array_out,1), ' does not match the data of size ',&
453  & (2*halo + 1 + var_in%iec - var_in%isc)
454  call mpp_error(fatal, trim(error_msg))
455  endif
456  i_off = 1 - (var_in%isc-halo)
457  endif
458 
459  ! Do error checking on the j-dimension and determine the array offsets.
460  if (present(jdim)) then
461  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
462  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
463  call mpp_error(fatal, trim(error_msg))
464  endif
465  if (size(array_out,2) /= (1+jdim(4)-jdim(1))) then
466  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
467  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_out,2)
468  call mpp_error(fatal, trim(error_msg))
469  endif
470  if ((var_in%jec-var_in%jsc) /= (jdim(3)-jdim(2)))&
471  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
472  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
473  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
474  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
475  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
476  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
477  & (2*halo + 1 + var_in%jec - var_in%jsc)
478  call mpp_error(fatal, trim(error_msg))
479  endif
480 
481  j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
482  else
483  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
484  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
485  & size(array_out,2), ' does not match the data of size ',&
486  & (2*halo + 1 + var_in%jec - var_in%jsc)
487  call mpp_error(fatal, trim(error_msg))
488  endif
489  j_off = 1 - (var_in%jsc-halo)
490  endif
491 
492  do j=var_in%jsc-halo,var_in%jec+halo
493  do i=var_in%isc-halo,var_in%iec+halo
494  array_out(i+i_off,j+j_off) = scale * var_in%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j)
495  enddo
496  enddo
497  end subroutine ct_extract_data_2d_
498 
499  !> @brief Extract a single k-level of a 3d field from a coupler_3d_bc_type
500  !!
501  !! Extract a single k-level of a 3-d field from a coupler_3d_bc_type into a two-dimensional array.
502  !!
503  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
504  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
505  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
506  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
507  !! @throw FATAL, "Disordered i-dimension index bound list"
508  !! @throw FATAL, "Disordered j-dimension index bound list"
509  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
510  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
511  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
512  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
513  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
514  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
515  !! @throw FATAL, "The extracted k-index of 'k' is outside of the valid range of 'ks' to 'ke'"
516  subroutine ct_extract_data_3d_2d_(var_in, bc_index, field_index, k_in, array_out,&
517  & scale_factor, halo_size, idim, jdim)
518  type(coupler_3d_bc_type), intent(in) :: var_in !< BC_type structure with the data to extract
519  integer, intent(in) :: bc_index !< The index of the boundary condition
520  !! that is being copied
521  integer, intent(in) :: field_index !< The index of the field in the
522  !! boundary condition that is being copied
523  integer, intent(in) :: k_in !< The k-index to extract
524  real(FMS_CP_KIND_), dimension(1:,1:), intent(out) :: array_out !< The recipient array for the field; its size
525  !! must match the size of the data being copied
526  !! unless idim and jdim are supplied.
527  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
528  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
529  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
530  !! the first dimension of the output array
531  !! in a non-decreasing list
532  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
533  !! the second dimension of the output array
534  !! in a non-decreasing list
535  character(len=*), parameter :: error_header =&
536  & '==>Error from coupler_types_mod (CT_extract_data_3d_2d):'
537  character(len=400) :: error_msg
538 
539  real(FMS_CP_KIND_) :: scale
540  integer :: i, j, halo, i_off, j_off
541  integer, parameter :: kindl = fms_cp_kind_
542 
543  if (bc_index <= 0) then
544  array_out(:,:) = 0.0_kindl
545  return
546  endif
547 
548  halo = 0
549  if (present(halo_size)) halo = halo_size
550  scale = 1.0_kindl
551  if (present(scale_factor)) scale = scale_factor
552 
553  if ((var_in%isc-var_in%isd < halo) .or. (var_in%ied-var_in%iec < halo))&
554  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
555  if ((var_in%jsc-var_in%jsd < halo) .or. (var_in%jed-var_in%jec < halo))&
556  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
557 
558  if (bc_index > var_in%num_bcs)&
559  & call mpp_error(fatal, trim(error_header)//" bc_index exceeds var_in%num_bcs.")
560  if (field_index > var_in%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
561  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
562  & trim(var_in%FMS_CP_BC_TYPE_(bc_index)%name) )
563 
564  ! Do error checking on the i-dimension and determine the array offsets.
565  if (present(idim)) then
566  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
567  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
568  call mpp_error(fatal, trim(error_msg))
569  endif
570  if (size(array_out,1) /= (1+idim(4)-idim(1))) then
571  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
572  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_out,1)
573  call mpp_error(fatal, trim(error_msg))
574  endif
575  if ((var_in%iec-var_in%isc) /= (idim(3)-idim(2)))&
576  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
577  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
578  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
579  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
580  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
581  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
582  & (2*halo + 1 + var_in%iec - var_in%isc)
583  call mpp_error(fatal, trim(error_msg))
584  endif
585 
586  i_off = (1-idim(1)) + (idim(2)-var_in%isc)
587  else
588  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
589  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
590  & size(array_out,1), ' does not match the data of size ',&
591  & (2*halo + 1 + var_in%iec - var_in%isc)
592  call mpp_error(fatal, trim(error_msg))
593  endif
594  i_off = 1 - (var_in%isc-halo)
595  endif
596 
597  ! Do error checking on the j-dimension and determine the array offsets.
598  if (present(jdim)) then
599  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
600  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
601  call mpp_error(fatal, trim(error_msg))
602  endif
603  if (size(array_out,2) /= (1+jdim(4)-jdim(1))) then
604  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
605  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_out,2)
606  call mpp_error(fatal, trim(error_msg))
607  endif
608  if ((var_in%jec-var_in%jsc) /= (jdim(3)-jdim(2)))&
609  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
610  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
611  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
612  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
613  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
614  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
615  & (2*halo + 1 + var_in%jec - var_in%jsc)
616  call mpp_error(fatal, trim(error_msg))
617  endif
618 
619  j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
620  else
621  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
622  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
623  & size(array_out,2), ' does not match the data of size ',&
624  & (2*halo + 1 + var_in%jec - var_in%jsc)
625  call mpp_error(fatal, trim(error_msg))
626  endif
627  j_off = 1 - (var_in%jsc-halo)
628  endif
629 
630  if ((k_in > var_in%ke) .or. (k_in < var_in%ks)) then
631  write (error_msg, *) trim(error_header), ' The extracted k-index of ', k_in,&
632  & ' is outside of the valid range of ', var_in%ks, ' to ', var_in%ke
633  call mpp_error(fatal, trim(error_msg))
634  endif
635 
636  do j=var_in%jsc-halo,var_in%jec+halo
637  do i=var_in%isc-halo,var_in%iec+halo
638  array_out(i+i_off,j+j_off) = scale * var_in%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j,k_in)
639  enddo
640  enddo
641  end subroutine ct_extract_data_3d_2d_
642 
643  !> @brief Extract single 3d field from a coupler_3d_bc_type
644  !!
645  !! Extract a single 3-d field from a coupler_3d_bc_type into a three-dimensional array.
646  !!
647  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
648  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
649  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
650  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
651  !! @throw FATAL, "Disordered i-dimension index bound list"
652  !! @throw FATAL, "Disordered j-dimension index bound list"
653  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
654  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
655  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
656  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
657  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
658  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
659  !! @throw FATAL, "The target array with k-dimension size 'n' does not match the data of size 'd'"
660  subroutine ct_extract_data_3d_(var_in, bc_index, field_index, array_out,&
661  & scale_factor, halo_size, idim, jdim)
662  type(coupler_3d_bc_type), intent(in) :: var_in !< BC_type structure with the data to extract
663  integer, intent(in) :: bc_index !< The index of the boundary condition
664  !! that is being copied
665  integer, intent(in) :: field_index !< The index of the field in the
666  !! boundary condition that is being copied
667  real(FMS_CP_KIND_), dimension(1:,1:,1:), intent(out) :: array_out !< The recipient array for the field; its size
668  !! must match the size of the data being copied
669  !! unless idim and jdim are supplied.
670  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
671  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
672  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
673  !! the first dimension of the output array
674  !! in a non-decreasing list
675  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
676  !! the second dimension of the output array
677  !! in a non-decreasing list
678 
679  character(len=*), parameter :: error_header =&
680  & '==>Error from coupler_types_mod (CT_extract_data_3d):'
681  character(len=400) :: error_msg
682 
683  real(FMS_CP_KIND_) :: scale
684  integer :: i, j, k, halo, i_off, j_off, k_off
685  integer, parameter :: kindl = fms_cp_kind_
686 
687  if (bc_index <= 0) then
688  array_out(:,:,:) = 0.0_kindl
689  return
690  endif
691 
692  halo = 0
693  if (present(halo_size)) halo = halo_size
694  scale = 1.0_kindl
695  if (present(scale_factor)) scale = scale_factor
696 
697  if ((var_in%isc-var_in%isd < halo) .or. (var_in%ied-var_in%iec < halo))&
698  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
699  if ((var_in%jsc-var_in%jsd < halo) .or. (var_in%jed-var_in%jec < halo))&
700  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
701 
702  if (bc_index > var_in%num_bcs)&
703  & call mpp_error(fatal, trim(error_header)//" bc_index exceeds var_in%num_bcs.")
704  if (field_index > var_in%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
705  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
706  & trim(var_in%FMS_CP_BC_TYPE_(bc_index)%name) )
707 
708  ! Do error checking on the i-dimension and determine the array offsets.
709  if (present(idim)) then
710  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
711  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
712  call mpp_error(fatal, trim(error_msg))
713  endif
714  if (size(array_out,1) /= (1+idim(4)-idim(1))) then
715  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
716  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_out,1)
717  call mpp_error(fatal, trim(error_msg))
718  endif
719  if ((var_in%iec-var_in%isc) /= (idim(3)-idim(2)))&
720  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
721  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
722  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
723  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
724  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
725  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
726  & (2*halo + 1 + var_in%iec - var_in%isc)
727  call mpp_error(fatal, trim(error_msg))
728  endif
729 
730  i_off = (1-idim(1)) + (idim(2)-var_in%isc)
731  else
732  if (size(array_out,1) < 2*halo + 1 + var_in%iec - var_in%isc) then
733  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
734  & size(array_out,1), ' does not match the data of size ',&
735  & (2*halo + 1 + var_in%iec - var_in%isc)
736  call mpp_error(fatal, trim(error_msg))
737  endif
738  i_off = 1 - (var_in%isc-halo)
739  endif
740 
741  ! Do error checking on the j-dimension and determine the array offsets.
742  if (present(jdim)) then
743  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
744  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
745  call mpp_error(fatal, trim(error_msg))
746  endif
747  if (size(array_out,2) /= (1+jdim(4)-jdim(1))) then
748  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
749  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_out,2)
750  call mpp_error(fatal, trim(error_msg))
751  endif
752  if ((var_in%jec-var_in%jsc) /= (jdim(3)-jdim(2)))&
753  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
754  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
755  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
756  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
757  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
758  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
759  & (2*halo + 1 + var_in%jec - var_in%jsc)
760  call mpp_error(fatal, trim(error_msg))
761  endif
762 
763  j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
764  else
765  if (size(array_out,2) < 2*halo + 1 + var_in%jec - var_in%jsc) then
766  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
767  & size(array_out,2), ' does not match the data of size ',&
768  & (2*halo + 1 + var_in%jec - var_in%jsc)
769  call mpp_error(fatal, trim(error_msg))
770  endif
771  j_off = 1 - (var_in%jsc-halo)
772  endif
773 
774  if (size(array_out,3) /= 1 + var_in%ke - var_in%ks) then
775  write (error_msg, *) trim(error_header), ' The target array with k-dimension size ',&
776  & size(array_out,3), ' does not match the data of size ',&
777  & (1 + var_in%ke - var_in%ks)
778  call mpp_error(fatal, trim(error_msg))
779  endif
780  k_off = 1 - var_in%ks
781 
782  do k=var_in%ks,var_in%ke
783  do j=var_in%jsc-halo,var_in%jec+halo
784  do i=var_in%isc-halo,var_in%iec+halo
785  array_out(i+i_off,j+j_off,k+k_off) = scale * var_in%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j,k)
786  enddo
787  enddo
788  enddo
789  end subroutine ct_extract_data_3d_
790 
791  !> @brief Set single 2d field in coupler_3d_bc_type
792  !!
793  !! Set a single 2-d field in a coupler_3d_bc_type from a two-dimensional array.
794  !!
795  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
796  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
797  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
798  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
799  !! @throw FATAL, "Disordered i-dimension index bound list"
800  !! @throw FATAL, "Disordered j-dimension index bound list"
801  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
802  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
803  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
804  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
805  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
806  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
807  subroutine ct_set_data_2d_(array_in, bc_index, field_index, var,&
808  & scale_factor, halo_size, idim, jdim)
809  real(FMS_CP_KIND_), dimension(1:,1:), intent(in) :: array_in !< The source array for the field; its size
810  !! must match the size of the data being copied
811  !! unless idim and jdim are supplied.
812  integer, intent(in) :: bc_index !< The index of the boundary condition
813  !! that is being copied
814  integer, intent(in) :: field_index !< The index of the field in the
815  !! boundary condition that is being copied
816  type(coupler_2d_bc_type), intent(inout) :: var !< BC_type structure with the data to set
817  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
818  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
819  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
820  !! the first dimension of the output array
821  !! in a non-decreasing list
822  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
823  !! the second dimension of the output array
824  !! in a non-decreasing list
825  character(len=*), parameter :: error_header =&
826  & '==>Error from coupler_types_mod (CT_set_data_2d):'
827  character(len=400) :: error_msg
828 
829  real(FMS_CP_KIND_) :: scale
830  integer :: i, j, halo, i_off, j_off
831  integer, parameter :: kindl = fms_cp_kind_
832 
833  if (bc_index <= 0) return
834 
835  halo = 0
836  if (present(halo_size)) halo = halo_size
837  scale = 1.0_kindl
838  if (present(scale_factor)) scale = scale_factor
839 
840  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
841  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
842  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
843  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
844 
845  if (bc_index > var%num_bcs) &
846  call mpp_error(fatal, trim(error_header)//" bc_index exceeds var%num_bcs.")
847  if (field_index > var%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
848  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
849  & trim(var%FMS_CP_BC_TYPE_(bc_index)%name) )
850 
851  ! Do error checking on the i-dimension and determine the array offsets.
852  if (present(idim)) then
853  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
854  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
855  call mpp_error(fatal, trim(error_msg))
856  endif
857  if (size(array_in,1) /= (1+idim(4)-idim(1))) then
858  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
859  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_in,1)
860  call mpp_error(fatal, trim(error_msg))
861  endif
862  if ((var%iec-var%isc) /= (idim(3)-idim(2)))&
863  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
864  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
865  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
866  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
867  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
868  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
869  & (2*halo + 1 + var%iec - var%isc)
870  call mpp_error(fatal, trim(error_msg))
871  endif
872 
873  i_off = (1-idim(1)) + (idim(2)-var%isc)
874  else
875  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
876  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
877  & size(array_in,1), ' does not match the data of size ',&
878  & (2*halo + 1 + var%iec - var%isc)
879  call mpp_error(fatal, trim(error_msg))
880  endif
881  i_off = 1 - (var%isc-halo)
882  endif
883 
884  ! Do error checking on the j-dimension and determine the array offsets.
885  if (present(jdim)) then
886  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
887  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
888  call mpp_error(fatal, trim(error_msg))
889  endif
890  if (size(array_in,2) /= (1+jdim(4)-jdim(1))) then
891  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
892  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_in,2)
893  call mpp_error(fatal, trim(error_msg))
894  endif
895  if ((var%jec-var%jsc) /= (jdim(3)-jdim(2)))&
896  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
897  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
898  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
899  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
900  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
901  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
902  & (2*halo + 1 + var%jec - var%jsc)
903  call mpp_error(fatal, trim(error_msg))
904  endif
905 
906  j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
907  else
908  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
909  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
910  & size(array_in,2), ' does not match the data of size ',&
911  & (2*halo + 1 + var%jec - var%jsc)
912  call mpp_error(fatal, trim(error_msg))
913  endif
914  j_off = 1 - (var%jsc-halo)
915  endif
916 
917  do j=var%jsc-halo,var%jec+halo
918  do i=var%isc-halo,var%iec+halo
919  var%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j) = scale * array_in(i+i_off,j+j_off)
920  enddo
921  enddo
922  end subroutine ct_set_data_2d_
923 
924  !> @brief Set one k-level of a single 3d field in a coupler_3d_bc_type
925  !!
926  !! This subroutine sets a one k-level of a single 3-d field in a coupler_3d_bc_type from a
927  !! two-dimensional array.
928  !!
929  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
930  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
931  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
932  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
933  !! @throw FATAL, "Disordered i-dimension index bound list"
934  !! @throw FATAL, "Disordered j-dimension index bound list"
935  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
936  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
937  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
938  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
939  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
940  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
941  !! @throw FATAL, "The k-index of 'k' is outside of the valid range of 'ks' to 'ke'"
942  subroutine ct_set_data_2d_3d_(array_in, bc_index, field_index, k_out, var,&
943  & scale_factor, halo_size, idim, jdim)
944  real(FMS_CP_KIND_), dimension(1:,1:), intent(in) :: array_in !< The source array for the field; its size
945  !! must match the size of the data being copied
946  !! unless idim and jdim are supplied.
947  integer, intent(in) :: bc_index !< The index of the boundary condition
948  !! that is being copied
949  integer, intent(in) :: field_index !< The index of the field in the
950  !! boundary condition that is being copied
951  integer, intent(in) :: k_out !< The k-index to set
952  type(coupler_3d_bc_type), intent(inout) :: var !< BC_type structure with the data to be set
953  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
954  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
955  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
956  !! the first dimension of the output array
957  !! in a non-decreasing list
958  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
959  !! the second dimension of the output array
960  !! in a non-decreasing list
961 
962  character(len=*), parameter :: error_header =&
963  & '==>Error from coupler_types_mod (CT_set_data_3d_2d):'
964  character(len=400) :: error_msg
965 
966  real(FMS_CP_KIND_) :: scale
967  integer :: i, j, halo, i_off, j_off
968  integer, parameter :: kindl = fms_cp_kind_
969 
970  if (bc_index <= 0) return
971 
972  halo = 0
973  if (present(halo_size)) halo = halo_size
974  scale = 1.0_kindl
975  if (present(scale_factor)) scale = scale_factor
976 
977  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
978  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
979  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
980  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
981 
982  if (bc_index > var%num_bcs)&
983  & call mpp_error(fatal, trim(error_header)//" bc_index exceeds var%num_bcs.")
984  if (field_index > var%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
985  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
986  & trim(var%FMS_CP_BC_TYPE_(bc_index)%name) )
987 
988  ! Do error checking on the i-dimension and determine the array offsets.
989  if (present(idim)) then
990  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
991  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
992  call mpp_error(fatal, trim(error_msg))
993  endif
994  if (size(array_in,1) /= (1+idim(4)-idim(1))) then
995  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
996  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_in,1)
997  call mpp_error(fatal, trim(error_msg))
998  endif
999  if ((var%iec-var%isc) /= (idim(3)-idim(2)))&
1000  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
1001  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
1002  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
1003  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
1004  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
1005  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
1006  & (2*halo + 1 + var%iec - var%isc)
1007  call mpp_error(fatal, trim(error_msg))
1008  endif
1009 
1010  i_off = (1-idim(1)) + (idim(2)-var%isc)
1011  else
1012  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
1013  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
1014  & size(array_in,1), ' does not match the data of size ',&
1015  & (2*halo + 1 + var%iec - var%isc)
1016  call mpp_error(fatal, trim(error_msg))
1017  endif
1018  i_off = 1 - (var%isc-halo)
1019  endif
1020 
1021  ! Do error checking on the j-dimension and determine the array offsets.
1022  if (present(jdim)) then
1023  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
1024  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
1025  call mpp_error(fatal, trim(error_msg))
1026  endif
1027  if (size(array_in,2) /= (1+jdim(4)-jdim(1))) then
1028  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
1029  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_in,2)
1030  call mpp_error(fatal, trim(error_msg))
1031  endif
1032  if ((var%jec-var%jsc) /= (jdim(3)-jdim(2)))&
1033  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
1034  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
1035  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
1036  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
1037  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
1038  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
1039  & (2*halo + 1 + var%jec - var%jsc)
1040  call mpp_error(fatal, trim(error_msg))
1041  endif
1042 
1043  j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
1044  else
1045  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
1046  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
1047  & size(array_in,2), ' does not match the data of size ',&
1048  & (2*halo + 1 + var%jec - var%jsc)
1049  call mpp_error(fatal, trim(error_msg))
1050  endif
1051  j_off = 1 - (var%jsc-halo)
1052  endif
1053 
1054  if ((k_out > var%ke) .or. (k_out < var%ks)) then
1055  write (error_msg, *) trim(error_header), ' The k-index of ', k_out,&
1056  & ' is outside of the valid range of ', var%ks, ' to ', var%ke
1057  call mpp_error(fatal, trim(error_msg))
1058  endif
1059 
1060  do j=var%jsc-halo,var%jec+halo
1061  do i=var%isc-halo,var%iec+halo
1062  var%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j,k_out) = scale * array_in(i+i_off,j+j_off)
1063  enddo
1064  enddo
1065  end subroutine ct_set_data_2d_3d_
1066 
1067  !> @brief Set a single 3d field in a coupler_3d_bc_type
1068  !!
1069  !! This subroutine sets a single 3-d field in a coupler_3d_bc_type from a three-dimensional array.
1070  !!
1071  !! @throw FATAL, "bc_index is present and exceeds var_in%num_bcs."
1072  !! @throw FATAL, "field_index exceeds num_fields for var_in%FMS_CP_BC_TYPE_(bc_incdx)%name"
1073  !! @throw FATAL, "Excessive i-direction halo size for the input structure."
1074  !! @throw FATAL, "Excessive j-direction halo size for the input structure."
1075  !! @throw FATAL, "Disordered i-dimension index bound list"
1076  !! @throw FATAL, "Disordered j-dimension index bound list"
1077  !! @throw FATAL, "The declared i-dimension size of 'n' does not match the actual size of 'a'"
1078  !! @throw FATAL, "The declared j-dimension size of 'n' does not match the actual size of 'a'"
1079  !! @throw FATAL, "There is an i-direction computational domain size mismatch."
1080  !! @throw FATAL, "There is an j-direction computational domain size mismatch."
1081  !! @throw FATAL, "The target array with i-dimension size 'n' is too small to match the data of size 'd'"
1082  !! @throw FATAL, "The target array with j-dimension size 'n' is too small to match the data of size 'd'"
1083  !! @throw FATAL, "The target array with K-dimension size 'n' is too small to match the data of size 'd'"
1084  subroutine ct_set_data_3d_(array_in, bc_index, field_index, var,&
1085  & scale_factor, halo_size, idim, jdim)
1086  real(FMS_CP_KIND_), dimension(1:,1:,1:), intent(in) :: array_in !< The source array for the field; its size
1087  !! must match the size of the data being copied
1088  !! unless idim and jdim are supplied.
1089  integer, intent(in) :: bc_index !< The index of the boundary condition
1090  !! that is being copied
1091  integer, intent(in) :: field_index !< The index of the field in the
1092  !! boundary condition that is being copied
1093  type(coupler_3d_bc_type), intent(inout) :: var !< BC_type structure with the data to be set
1094  real(FMS_CP_KIND_), optional, intent(in) :: scale_factor !< A scaling factor for the data that is being added
1095  integer, optional, intent(in) :: halo_size !< The extent of the halo to copy; 0 by default
1096  integer, dimension(4), optional, intent(in) :: idim !< The data and computational domain extents of
1097  !! the first dimension of the output array
1098  !! in a non-decreasing list
1099  integer, dimension(4), optional, intent(in) :: jdim !< The data and computational domain extents of
1100  !! the second dimension of the output array
1101  !! in a non-decreasing list
1102 
1103  character(len=*), parameter :: error_header =&
1104  & '==>Error from coupler_types_mod (CT_set_data_3d):'
1105  character(len=400) :: error_msg
1106 
1107  real(FMS_CP_KIND_) :: scale
1108  integer :: i, j, k, halo, i_off, j_off, k_off
1109  integer, parameter :: kindl = fms_cp_kind_
1110 
1111  if (bc_index <= 0) return
1112 
1113  halo = 0
1114  if (present(halo_size)) halo = halo_size
1115  scale = 1.0_kindl
1116  if (present(scale_factor)) scale = scale_factor
1117 
1118  if ((var%isc-var%isd < halo) .or. (var%ied-var%iec < halo))&
1119  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the input structure.")
1120  if ((var%jsc-var%jsd < halo) .or. (var%jed-var%jec < halo))&
1121  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the input structure.")
1122 
1123  if (bc_index > var%num_bcs)&
1124  & call mpp_error(fatal, trim(error_header)//" bc_index exceeds var%num_bcs.")
1125  if (field_index > var%FMS_CP_BC_TYPE_(bc_index)%num_fields)&
1126  & call mpp_error(fatal, trim(error_header)//" field_index exceeds num_fields for" //&
1127  & trim(var%FMS_CP_BC_TYPE_(bc_index)%name) )
1128 
1129  ! Do error checking on the i-dimension and determine the array offsets.
1130  if (present(idim)) then
1131  if ((idim(1) > idim(2)) .or. (idim(3) > idim(4))) then
1132  write (error_msg, *) trim(error_header), ' Disordered i-dimension index bound list ', idim
1133  call mpp_error(fatal, trim(error_msg))
1134  endif
1135  if (size(array_in,1) /= (1+idim(4)-idim(1))) then
1136  write (error_msg, *) trim(error_header), ' The declared i-dimension size of ',&
1137  & (1+idim(4)-idim(1)), ' does not match the actual size of ', size(array_in,1)
1138  call mpp_error(fatal, trim(error_msg))
1139  endif
1140  if ((var%iec-var%isc) /= (idim(3)-idim(2)))&
1141  & call mpp_error(fatal, trim(error_header)//" There is an i-direction computational domain size mismatch.")
1142  if ((idim(2)-idim(1) < halo) .or. (idim(4)-idim(3) < halo))&
1143  & call mpp_error(fatal, trim(error_header)//" Excessive i-direction halo size for the output array.")
1144  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
1145  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
1146  & (1+idim(4)-idim(1)), ' is too small to match the data of size ',&
1147  & (2*halo + 1 + var%iec - var%isc)
1148  call mpp_error(fatal, trim(error_msg))
1149  endif
1150 
1151  i_off = (1-idim(1)) + (idim(2)-var%isc)
1152  else
1153  if (size(array_in,1) < 2*halo + 1 + var%iec - var%isc) then
1154  write (error_msg, *) trim(error_header), ' The target array with i-dimension size ',&
1155  & size(array_in,1), ' does not match the data of size ',&
1156  & (2*halo + 1 + var%iec - var%isc)
1157  call mpp_error(fatal, trim(error_msg))
1158  endif
1159  i_off = 1 - (var%isc-halo)
1160  endif
1161 
1162  ! Do error checking on the j-dimension and determine the array offsets.
1163  if (present(jdim)) then
1164  if ((jdim(1) > jdim(2)) .or. (jdim(3) > jdim(4))) then
1165  write (error_msg, *) trim(error_header), ' Disordered j-dimension index bound list ', jdim
1166  call mpp_error(fatal, trim(error_msg))
1167  endif
1168  if (size(array_in,2) /= (1+jdim(4)-jdim(1))) then
1169  write (error_msg, *) trim(error_header), ' The declared j-dimension size of ',&
1170  & (1+jdim(4)-jdim(1)), ' does not match the actual size of ', size(array_in,2)
1171  call mpp_error(fatal, trim(error_msg))
1172  endif
1173  if ((var%jec-var%jsc) /= (jdim(3)-jdim(2)))&
1174  & call mpp_error(fatal, trim(error_header)//" There is an j-direction computational domain size mismatch.")
1175  if ((jdim(2)-jdim(1) < halo) .or. (jdim(4)-jdim(3) < halo))&
1176  & call mpp_error(fatal, trim(error_header)//" Excessive j-direction halo size for the output array.")
1177  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
1178  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
1179  & (1+jdim(4)-jdim(1)), ' is too small to match the data of size ',&
1180  & (2*halo + 1 + var%jec - var%jsc)
1181  call mpp_error(fatal, trim(error_msg))
1182  endif
1183 
1184  j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
1185  else
1186  if (size(array_in,2) < 2*halo + 1 + var%jec - var%jsc) then
1187  write (error_msg, *) trim(error_header), ' The target array with j-dimension size ',&
1188  & size(array_in,2), ' does not match the data of size ',&
1189  & (2*halo + 1 + var%jec - var%jsc)
1190  call mpp_error(fatal, trim(error_msg))
1191  endif
1192  j_off = 1 - (var%jsc-halo)
1193  endif
1194 
1195  if (size(array_in,3) /= 1 + var%ke - var%ks) then
1196  write (error_msg, *) trim(error_header), ' The target array with k-dimension size ',&
1197  & size(array_in,3), ' does not match the data of size ',&
1198  & (1 + var%ke - var%ks)
1199  call mpp_error(fatal, trim(error_msg))
1200  endif
1201  k_off = 1 - var%ks
1202 
1203  do k=var%ks,var%ke
1204  do j=var%jsc-halo,var%jec+halo
1205  do i=var%isc-halo,var%iec+halo
1206  var%FMS_CP_BC_TYPE_(bc_index)%field(field_index)%values(i,j,k) = scale * array_in(i+i_off,j+j_off,k+k_off)
1207  enddo
1208  enddo
1209  enddo
1210  end subroutine ct_set_data_3d_
subroutine ct_extract_data_2d_(var_in, bc_index, field_index, array_out, scale_factor, halo_size, idim, jdim)
Extract a 2d field from a coupler_2d_bc_type.
subroutine ct_rescale_data_2d_(var, scale, halo_size, bc_index, field_index, exclude_flux_type, only_flux_type, pass_through_ice)
Rescales the fields in the fields in the elements of a coupler_2d_bc_type.
subroutine ct_set_data_3d_(array_in, bc_index, field_index, var, scale_factor, halo_size, idim, jdim)
Set a single 3d field in a coupler_3d_bc_type.
subroutine ct_extract_data_3d_2d_(var_in, bc_index, field_index, k_in, array_out, scale_factor, halo_size, idim, jdim)
Extract a single k-level of a 3d field from a coupler_3d_bc_type.
subroutine ct_set_data_2d_3d_(array_in, bc_index, field_index, k_out, var, scale_factor, halo_size, idim, jdim)
Set one k-level of a single 3d field in a coupler_3d_bc_type.
subroutine ct_rescale_data_3d_(var, scale, halo_size, bc_index, field_index, exclude_flux_type, only_flux_type, pass_through_ice)
Rescales the fields in the elements of a coupler_3d_bc_type.
subroutine ct_increment_data_2d_3d_(var_in, weights, var, halo_size, bc_index, field_index, scale_factor, scale_prev, exclude_flux_type, only_flux_type, pass_through_ice)
Increment data in the elements of a coupler_2d_bc_type with weighted averages of elements of a couple...
subroutine ct_extract_data_3d_(var_in, bc_index, field_index, array_out, scale_factor, halo_size, idim, jdim)
Extract single 3d field from a coupler_3d_bc_type.
subroutine ct_set_data_2d_(array_in, bc_index, field_index, var, scale_factor, halo_size, idim, jdim)
Set single 2d field in coupler_3d_bc_type.