35 & exclude_flux_type, only_flux_type, pass_through_ice)
36 type(coupler_2d_bc_type),
intent(inout) :: var
37 real(FMS_CP_KIND_),
intent(in) :: scale
38 integer,
optional,
intent(in) :: halo_size
40 integer,
optional,
intent(in) :: bc_index
42 integer,
optional,
intent(in) :: field_index
44 character(len=*),
optional,
intent(in) :: exclude_flux_type
46 character(len=*),
optional,
intent(in) :: only_flux_type
48 logical,
optional,
intent(in) :: pass_through_ice
52 integer :: i, j, m, n, n1, n2, halo
53 integer,
parameter :: kindl = fms_cp_kind_
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) )
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.")
67 if (
present(halo_size)) halo = halo_size
71 if (
present(bc_index))
then
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.")
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)
94 do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
95 if (
present(field_index))
then
96 if (m /= field_index) cycle
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
107 var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:) = 0.0_kindl
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)
126 & exclude_flux_type, only_flux_type, pass_through_ice)
127 type(coupler_3d_bc_type),
intent(inout) :: var
128 real(FMS_CP_KIND_),
intent(in) :: scale
129 integer,
optional,
intent(in) :: halo_size
131 integer,
optional,
intent(in) :: bc_index
133 integer,
optional,
intent(in) :: field_index
135 character(len=*),
optional,
intent(in) :: exclude_flux_type
137 character(len=*),
optional,
intent(in) :: only_flux_type
139 logical,
optional,
intent(in) :: pass_through_ice
143 integer :: i, j, k, m, n, n1, n2, halo
144 integer,
parameter :: kindl = fms_cp_kind_
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) )
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.")
158 if (
present(halo_size)) halo = halo_size
162 if (
present(bc_index))
then
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.")
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
185 do m = 1, var%FMS_CP_BC_TYPE_(n)%num_fields
186 if (
present(field_index))
then
187 if (m /= field_index) cycle
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
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
200 var%FMS_CP_BC_TYPE_(n)%field(m)%values(:,:,:) = 0.0_kindl
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)
234 & scale_factor, scale_prev, exclude_flux_type, only_flux_type, pass_through_ice)
235 type(coupler_3d_bc_type),
intent(in) :: var_in
236 real(FMS_CP_KIND_),
dimension(:,:,:),
intent(in) :: weights
240 type(coupler_2d_bc_type),
intent(inout) :: var
241 integer,
optional,
intent(in) :: halo_size
242 integer,
optional,
intent(in) :: bc_index
244 integer,
optional,
intent(in) :: field_index
246 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
248 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_prev
250 character(len=*),
optional,
intent(in) :: exclude_flux_type
252 character(len=*),
optional,
intent(in) :: only_flux_type
254 logical,
optional,
intent(in) :: pass_through_ice
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
261 integer,
parameter :: kindl = fms_cp_kind_
264 if (
present(scale_factor)) scale = scale_factor
266 if (
present(scale_prev)) sc_prev = scale_prev
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) )
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.")
280 if (
present(halo_size)) halo = halo_size
284 if (
present(bc_index))
then
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.")
309 if ((1+var%iec-var%isc) ==
size(weights,1))
then
311 elseif ((1+var%ied-var%isd) ==
size(weights,1))
then
313 elseif ((1+var_in%ied-var_in%isd) ==
size(weights,1))
then
314 iow = 1 + (var_in%isc - var_in%isd) - var%isc
316 call mpp_error(fatal, &
317 &
"CT_increment_data_2d_3d: weights array must be the i-size of a computational or data domain.")
319 if ((1+var%jec-var%jsc) ==
size(weights,2))
then
321 elseif ((1+var%jed-var%jsd) ==
size(weights,2))
then
323 elseif ((1+var_in%jed-var_in%jsd) ==
size(weights,2))
then
324 jow = 1 + (var_in%jsc - var_in%jsd) - var%jsc
326 call mpp_error(fatal, &
327 &
"CT_increment_data_2d_3d: weights array must be the j-size of a computational or data domain.")
330 io1 = var_in%isc - var%isc
331 jo1 = var_in%jsc - var%jsc
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
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
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)
380 & scale_factor, halo_size, idim, jdim)
381 type(coupler_2d_bc_type),
intent(in) :: var_in
382 integer,
intent(in) :: bc_index
384 integer,
intent(in) :: field_index
386 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(out) :: array_out
389 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
390 integer,
optional,
intent(in) :: halo_size
391 integer,
dimension(4),
optional,
intent(in) :: idim
394 integer,
dimension(4),
optional,
intent(in) :: jdim
398 character(len=*),
parameter :: error_header =&
399 &
'==>Error from coupler_types_mod (CT_extract_data_2d):'
400 character(len=400) :: error_msg
402 real(FMS_CP_KIND_) :: scale
403 integer :: i, j, halo, i_off, j_off
404 integer,
parameter :: kindl = fms_cp_kind_
406 if (bc_index <= 0)
then
407 array_out(:,:) = 0.0_kindl
412 if (
present(halo_size)) halo = halo_size
414 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
449 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
457 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
482 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
490 j_off = 1 - (var_in%jsc-halo)
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)
518 & scale_factor, halo_size, idim, jdim)
519 type(coupler_3d_bc_type),
intent(in) :: var_in
520 integer,
intent(in) :: bc_index
522 integer,
intent(in) :: field_index
524 integer,
intent(in) :: k_in
525 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(out) :: array_out
528 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
529 integer,
optional,
intent(in) :: halo_size
530 integer,
dimension(4),
optional,
intent(in) :: idim
533 integer,
dimension(4),
optional,
intent(in) :: jdim
536 character(len=*),
parameter :: error_header =&
537 &
'==>Error from coupler_types_mod (CT_extract_data_3d_2d):'
538 character(len=400) :: error_msg
540 real(FMS_CP_KIND_) :: scale
541 integer :: i, j, halo, i_off, j_off
542 integer,
parameter :: kindl = fms_cp_kind_
544 if (bc_index <= 0)
then
545 array_out(:,:) = 0.0_kindl
550 if (
present(halo_size)) halo = halo_size
552 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
587 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
595 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
620 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
628 j_off = 1 - (var_in%jsc-halo)
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))
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)
662 & scale_factor, halo_size, idim, jdim)
663 type(coupler_3d_bc_type),
intent(in) :: var_in
664 integer,
intent(in) :: bc_index
666 integer,
intent(in) :: field_index
668 real(FMS_CP_KIND_),
dimension(1:,1:,1:),
intent(out) :: array_out
671 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
672 integer,
optional,
intent(in) :: halo_size
673 integer,
dimension(4),
optional,
intent(in) :: idim
676 integer,
dimension(4),
optional,
intent(in) :: jdim
680 character(len=*),
parameter :: error_header =&
681 &
'==>Error from coupler_types_mod (CT_extract_data_3d):'
682 character(len=400) :: error_msg
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_
688 if (bc_index <= 0)
then
689 array_out(:,:,:) = 0.0_kindl
694 if (
present(halo_size)) halo = halo_size
696 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
731 i_off = (1-idim(1)) + (idim(2)-var_in%isc)
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))
739 i_off = 1 - (var_in%isc-halo)
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))
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))
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))
764 j_off = (1-jdim(1)) + (jdim(2)-var_in%jsc)
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))
772 j_off = 1 - (var_in%jsc-halo)
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))
781 k_off = 1 - var_in%ks
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)
809 & scale_factor, halo_size, idim, jdim)
810 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(in) :: array_in
813 integer,
intent(in) :: bc_index
815 integer,
intent(in) :: field_index
817 type(coupler_2d_bc_type),
intent(inout) :: var
818 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
819 integer,
optional,
intent(in) :: halo_size
820 integer,
dimension(4),
optional,
intent(in) :: idim
823 integer,
dimension(4),
optional,
intent(in) :: jdim
826 character(len=*),
parameter :: error_header =&
827 &
'==>Error from coupler_types_mod (CT_set_data_2d):'
828 character(len=400) :: error_msg
830 real(FMS_CP_KIND_) :: scale
831 integer :: i, j, halo, i_off, j_off
832 integer,
parameter :: kindl = fms_cp_kind_
834 if (bc_index <= 0)
return
837 if (
present(halo_size)) halo = halo_size
839 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
874 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
882 i_off = 1 - (var%isc-halo)
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))
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))
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))
907 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
915 j_off = 1 - (var%jsc-halo)
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)
944 & scale_factor, halo_size, idim, jdim)
945 real(FMS_CP_KIND_),
dimension(1:,1:),
intent(in) :: array_in
948 integer,
intent(in) :: bc_index
950 integer,
intent(in) :: field_index
952 integer,
intent(in) :: k_out
953 type(coupler_3d_bc_type),
intent(inout) :: var
954 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
955 integer,
optional,
intent(in) :: halo_size
956 integer,
dimension(4),
optional,
intent(in) :: idim
959 integer,
dimension(4),
optional,
intent(in) :: jdim
963 character(len=*),
parameter :: error_header =&
964 &
'==>Error from coupler_types_mod (CT_set_data_3d_2d):'
965 character(len=400) :: error_msg
967 real(FMS_CP_KIND_) :: scale
968 integer :: i, j, halo, i_off, j_off
969 integer,
parameter :: kindl = fms_cp_kind_
971 if (bc_index <= 0)
return
974 if (
present(halo_size)) halo = halo_size
976 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
1011 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
1019 i_off = 1 - (var%isc-halo)
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))
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))
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))
1044 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
1052 j_off = 1 - (var%jsc-halo)
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))
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)
1086 & scale_factor, halo_size, idim, jdim)
1087 real(FMS_CP_KIND_),
dimension(1:,1:,1:),
intent(in) :: array_in
1090 integer,
intent(in) :: bc_index
1092 integer,
intent(in) :: field_index
1094 type(coupler_3d_bc_type),
intent(inout) :: var
1095 real(FMS_CP_KIND_),
optional,
intent(in) :: scale_factor
1096 integer,
optional,
intent(in) :: halo_size
1097 integer,
dimension(4),
optional,
intent(in) :: idim
1100 integer,
dimension(4),
optional,
intent(in) :: jdim
1104 character(len=*),
parameter :: error_header =&
1105 &
'==>Error from coupler_types_mod (CT_set_data_3d):'
1106 character(len=400) :: error_msg
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_
1112 if (bc_index <= 0)
return
1115 if (
present(halo_size)) halo = halo_size
1117 if (
present(scale_factor)) scale = scale_factor
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.")
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) )
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))
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))
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))
1152 i_off = (1-idim(1)) + (idim(2)-var%isc)
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))
1160 i_off = 1 - (var%isc-halo)
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))
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))
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))
1185 j_off = (1-jdim(1)) + (jdim(2)-var%jsc)
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))
1193 j_off = 1 - (var%jsc-halo)
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))
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)
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.