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