34 & exclude_flux_type, only_flux_type, pass_through_ice)
35 type(coupler_2d_bc_type),
intent(inout) :: var
36 real(FMS_CP_KIND_),
intent(in) :: scale
37 integer,
optional,
intent(in) :: halo_size
39 integer,
optional,
intent(in) :: bc_index
41 integer,
optional,
intent(in) :: field_index
43 character(len=*),
optional,
intent(in) :: exclude_flux_type
45 character(len=*),
optional,
intent(in) :: only_flux_type
47 logical,
optional,
intent(in) :: pass_through_ice
51 integer :: i, j, m, n, n1, n2, halo
52 integer,
parameter :: kindl = fms_cp_kind_
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) )
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.")
66 if (
present(halo_size)) halo = halo_size
70 if (
present(bc_index))
then
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.")
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)
93 do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
94 if (
present(field_index))
then
95 if (m /= field_index) cycle
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
106 var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:) = 0.0_kindl
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)
125 & exclude_flux_type, only_flux_type, pass_through_ice)
126 type(coupler_3d_bc_type),
intent(inout) :: var
127 real(FMS_CP_KIND_),
intent(in) :: scale
128 integer,
optional,
intent(in) :: halo_size
130 integer,
optional,
intent(in) :: bc_index
132 integer,
optional,
intent(in) :: field_index
134 character(len=*),
optional,
intent(in) :: exclude_flux_type
136 character(len=*),
optional,
intent(in) :: only_flux_type
138 logical,
optional,
intent(in) :: pass_through_ice
142 integer :: i, j, k, m, n, n1, n2, halo
143 integer,
parameter :: kindl = fms_cp_kind_
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) )
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.")
157 if (
present(halo_size)) halo = halo_size
161 if (
present(bc_index))
then
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.")
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
184 do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
185 if (
present(field_index))
then
186 if (m /= field_index) cycle
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
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
199 var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:,:) = 0.0_kindl
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)
233 & scale_factor, scale_prev, exclude_flux_type, only_flux_type, pass_through_ice)
234 type(coupler_3d_bc_type),
intent(in) :: var_in
235 real(FMS_CP_KIND_),
dimension(:,:,:),
intent(in) :: weights
239 type(coupler_2d_bc_type),
intent(inout) :: var
240 integer,
optional,
intent(in) :: halo_size
241 integer,
optional,
intent(in) :: bc_index
243 integer,
optional,
intent(in) :: field_index
245 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
247 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_prev
249 character(len=*),
optional,
intent(in) :: exclude_flux_type
251 character(len=*),
optional,
intent(in) :: only_flux_type
253 logical,
optional,
intent(in) :: pass_through_ice
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
260 integer,
parameter :: kindl = fms_cp_kind_
263 if (
present(scale_factor)) scale = scale_factor
265 if (
present(scale_prev)) sc_prev = scale_prev
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) )
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.")
279 if (
present(halo_size)) halo = halo_size
283 if (
present(bc_index))
then
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.")
308 if ((1+var%iec-var%isc) ==
size(weights,1))
then
310 elseif ((1+var%ied-var%isd) ==
size(weights,1))
then
312 elseif ((1+var_in%ied-var_in%isd) ==
size(weights,1))
then
313 iow = 1 + (var_in%isc - var_in%isd) - var%isc
315 call mpp_error(fatal, &
316 &
"CT_increment_data_2d_3d: weights array must be the i-size of a computational or data domain.")
318 if ((1+var%jec-var%jsc) ==
size(weights,2))
then
320 elseif ((1+var%jed-var%jsd) ==
size(weights,2))
then
322 elseif ((1+var_in%jed-var_in%jsd) ==
size(weights,2))
then
323 jow = 1 + (var_in%jsc - var_in%jsd) - var%jsc
325 call mpp_error(fatal, &
326 &
"CT_increment_data_2d_3d: weights array must be the j-size of a computational or data domain.")
329 io1 = var_in%isc - var%isc
330 jo1 = var_in%jsc - var%jsc
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
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
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)
379 & scale_factor, halo_size, idim, jdim)
380 type(coupler_2d_bc_type),
intent(in) :: var_in
381 integer,
intent(in) :: bc_index
383 integer,
intent(in) :: field_index
385 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(out) :: array_out
388 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
389 integer,
optional,
intent(in) :: halo_size
390 integer,
dimension(4),
optional,
intent(in) :: idim
393 integer,
dimension(4),
optional,
intent(in) :: jdim
397 character(len=*),
parameter :: error_header =&
398 &
'==>Error from coupler_types_mod (CT_extract_data_2d):'
399 character(len=400) :: error_msg
401 real(FMS_CP_KIND_) :: scale
402 integer :: i, j, halo, i_off, j_off
403 integer,
parameter :: kindl = fms_cp_kind_
405 if (bc_index <= 0)
then
406 array_out(:,:) = 0.0_kindl
411 if (
present(halo_size)) halo = halo_size
413 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
448 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
456 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
481 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
489 j_off = 1 - (var_in%jsc-halo)
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)
517 & scale_factor, halo_size, idim, jdim)
518 type(coupler_3d_bc_type),
intent(in) :: var_in
519 integer,
intent(in) :: bc_index
521 integer,
intent(in) :: field_index
523 integer,
intent(in) :: k_in
524 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(out) :: array_out
527 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
528 integer,
optional,
intent(in) :: halo_size
529 integer,
dimension(4),
optional,
intent(in) :: idim
532 integer,
dimension(4),
optional,
intent(in) :: jdim
535 character(len=*),
parameter :: error_header =&
536 &
'==>Error from coupler_types_mod (CT_extract_data_3d_2d):'
537 character(len=400) :: error_msg
539 real(FMS_CP_KIND_) :: scale
540 integer :: i, j, halo, i_off, j_off
541 integer,
parameter :: kindl = fms_cp_kind_
543 if (bc_index <= 0)
then
544 array_out(:,:) = 0.0_kindl
549 if (
present(halo_size)) halo = halo_size
551 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
586 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
594 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
619 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
627 j_off = 1 - (var_in%jsc-halo)
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))
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)
661 & scale_factor, halo_size, idim, jdim)
662 type(coupler_3d_bc_type),
intent(in) :: var_in
663 integer,
intent(in) :: bc_index
665 integer,
intent(in) :: field_index
667 real(FMS_CP_KIND_),
dimension(1:,1:,1:),
intent(out) :: array_out
670 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
671 integer,
optional,
intent(in) :: halo_size
672 integer,
dimension(4),
optional,
intent(in) :: idim
675 integer,
dimension(4),
optional,
intent(in) :: jdim
679 character(len=*),
parameter :: error_header =&
680 &
'==>Error from coupler_types_mod (CT_extract_data_3d):'
681 character(len=400) :: error_msg
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_
687 if (bc_index <= 0)
then
688 array_out(:,:,:) = 0.0_kindl
693 if (
present(halo_size)) halo = halo_size
695 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
730 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
738 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
763 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
771 j_off = 1 - (var_in%jsc-halo)
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))
780 k_off = 1 - var_in%ks
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)
808 & scale_factor, halo_size, idim, jdim)
809 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(in) :: array_in
812 integer,
intent(in) :: bc_index
814 integer,
intent(in) :: field_index
816 type(coupler_2d_bc_type),
intent(inout) :: var
817 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
818 integer,
optional,
intent(in) :: halo_size
819 integer,
dimension(4),
optional,
intent(in) :: idim
822 integer,
dimension(4),
optional,
intent(in) :: jdim
825 character(len=*),
parameter :: error_header =&
826 &
'==>Error from coupler_types_mod (CT_set_data_2d):'
827 character(len=400) :: error_msg
829 real(FMS_CP_KIND_) :: scale
830 integer :: i, j, halo, i_off, j_off
831 integer,
parameter :: kindl = fms_cp_kind_
833 if (bc_index <= 0)
return
836 if (
present(halo_size)) halo = halo_size
838 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
873 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
881 i_off = 1 - (var%isc-halo)
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))
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))
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))
906 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
914 j_off = 1 - (var%jsc-halo)
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)
943 & scale_factor, halo_size, idim, jdim)
944 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(in) :: array_in
947 integer,
intent(in) :: bc_index
949 integer,
intent(in) :: field_index
951 integer,
intent(in) :: k_out
952 type(coupler_3d_bc_type),
intent(inout) :: var
953 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
954 integer,
optional,
intent(in) :: halo_size
955 integer,
dimension(4),
optional,
intent(in) :: idim
958 integer,
dimension(4),
optional,
intent(in) :: jdim
962 character(len=*),
parameter :: error_header =&
963 &
'==>Error from coupler_types_mod (CT_set_data_3d_2d):'
964 character(len=400) :: error_msg
966 real(FMS_CP_KIND_) :: scale
967 integer :: i, j, halo, i_off, j_off
968 integer,
parameter :: kindl = fms_cp_kind_
970 if (bc_index <= 0)
return
973 if (
present(halo_size)) halo = halo_size
975 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
1010 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
1018 i_off = 1 - (var%isc-halo)
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))
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))
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))
1043 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
1051 j_off = 1 - (var%jsc-halo)
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))
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)
1085 & scale_factor, halo_size, idim, jdim)
1086 real(FMS_CP_KIND_),
dimension(1:,1:,1:),
intent(in) :: array_in
1089 integer,
intent(in) :: bc_index
1091 integer,
intent(in) :: field_index
1093 type(coupler_3d_bc_type),
intent(inout) :: var
1094 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
1095 integer,
optional,
intent(in) :: halo_size
1096 integer,
dimension(4),
optional,
intent(in) :: idim
1099 integer,
dimension(4),
optional,
intent(in) :: jdim
1103 character(len=*),
parameter :: error_header =&
1104 &
'==>Error from coupler_types_mod (CT_set_data_3d):'
1105 character(len=400) :: error_msg
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_
1111 if (bc_index <= 0)
return
1114 if (
present(halo_size)) halo = halo_size
1116 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
1151 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
1159 i_off = 1 - (var%isc-halo)
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))
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))
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))
1184 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
1192 j_off = 1 - (var%jsc-halo)
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))
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)
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.