psmile_merge_fields.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2010, DKRZ, Hamburg, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !
00006 ! !DESCRIPTION:
00007 !
00008 ! The routine psmile_merge_fields postprocesses the results of enddef
00009 ! process. It collects the send and receive information generated by
00010 ! the search process of all block for transient and merges them into a
00011 ! single field, which is easier to handle than having a field for each
00012 ! block.
00013 !  In order to do this the internal data structures need to be edited.
00014 ! This is an evil hack, VERY error-prone, and VERY hard to maintain.
00015 ! Hence, I do not like it at all. But I currently see no better way.
00016 !  After the merge the user_var_ids generated by psmile_store_user_data_var
00017 ! match the internal psmile_var_ids. Therefore, no translation is needed
00018 ! (e.g. in psmile_put or get)
00019 !
00020 !
00021 ! !REVISION HISTORY:
00022 !
00023 !   Date      Programmer   Description
00024 ! ----------  ----------   -----------
00025 ! 08.11.10    M. Hanke     created
00026 !
00027 !----------------------------------------------------------------------
00028 !
00029 !  $Id: psmile_merge_fields.F90 3069 2011-03-25 14:01:22Z redler $
00030 !  $Author: redler $
00031 !
00032 !----------------------------------------------------------------------
00033 
00034 subroutine psmile_merge_fields (ierror)
00035 
00036    use psmile, only            : GridFunction, Method, Mask, Fields, Methods, Masks, &
00037                                  Grids, taskout_type, taskin_type,                   &
00038                                  send_field_information, recv_field_information
00039    use psmile_user_data, only  : var_data, var_id_map,     &
00040                                  point_data, point_id_map, &
00041                                  grid_data, grid_id_map,   &
00042                                  mask_data, mask_id_map
00043    use psmile_grid, only       : get_size_of_shape
00044    use psmile_multimap, only   : get_num_values, get_values, get_value, &
00045                                  multimap, init_multimap, add_pair,     &
00046                                  is_valid_key, free_multimap
00047    use psmile_common, only     : number_of_fields_allocated, psmile_status_free, &
00048                                  Number_of_Methods_allocated,                    &
00049                                  Number_of_Fields_allocated,                     &
00050                                  Number_of_Masks_allocated,                      &
00051                                  ch_id, psmile_undef, ndim_3d
00052    use psmile_reallocate, only : psmile_realloc
00053    use prism_constants, only   : prism_undefined
00054 
00055    implicit none
00056 
00057    integer, intent(out) :: ierror
00058 
00059    integer                     :: i, j, k
00060    integer                     :: num_blocks
00061    integer, allocatable        :: var_ids(:), method_ids(:), mask_ids(:)
00062    type(GridFunction), pointer :: field_final(:)
00063    type(Method), pointer       :: method_final(:)
00064    type(Mask), pointer         :: mask_final(:)
00065 
00066    type info_index_offset
00067       integer :: send_direct, send_coupler, send_appl, recv_direct, recv_coupler
00068    end type
00069 
00070    type method_patch_info
00071       type (info_index_offset) :: info_offset
00072       integer :: new_method_id
00073       integer :: new_mask_id
00074    end type
00075 
00076    type (info_index_offset)              :: info_offset
00077    type (method_patch_info), allocatable :: method_patches (:)
00078    type (multimap)                       :: old_to_new_var_id
00079 
00080 #ifdef VERBOSE
00081    print 9990, trim(ch_id)
00082    call psmile_flushstd
00083 #endif /* VERBOSE */
00084 
00085    ierror = 0
00086 
00087    call init_multimap (old_to_new_var_id)
00088 
00089    if ( associated (grid_data)) &
00090       allocate (var_ids   (get_max_num_blocks ()), &
00091                 method_ids(get_max_num_blocks ()), &
00092                 mask_ids  (get_max_num_blocks ()))
00093 
00094    if (number_of_fields_allocated > 0) &
00095       allocate (method_patches(number_of_fields_allocated))
00096 
00097 !===============================================================================
00098 !  1st merge methods/points
00099 !===============================================================================
00100 
00101    if (associated (point_data)) then
00102 
00103       allocate (method_final(size (point_data)))
00104 
00105       ! for all point sets
00106       do i = 1, size (point_data)
00107 
00108          num_blocks = get_num_values (point_id_map, i)
00109 
00110          ! get the method ids of the methods associated of the current method id
00111          method_ids (1:num_blocks) = get_values (point_id_map, i, num_blocks)
00112 
00113          do j = 1, num_blocks
00114 
00115             if (j == 1) then
00116 
00117                ! make a copy of the first method
00118                call copy_method (method_final(i), methods(method_ids(1)))
00119                info_offset = info_index_offset(0,0,0,0,0)
00120 
00121             else
00122 
00123                ! merge the current method with the first one
00124                call merge_method (method_final(i), methods(method_ids(j)), info_offset)
00125             endif
00126 
00127             ! for all fields
00128             do k = 1, number_of_fields_allocated
00129 
00130                ! find fields that reference the current method
00131                if (Fields(k)%status /= psmile_status_free .and. &
00132                    Fields(k)%method_id == method_ids(j)) then
00133 
00134                   ! this information is later neede to adjust the merged fields to the new method structure
00135                   method_patches(k)%info_offset = info_offset
00136                   method_patches(k)%new_method_id = i
00137 
00138                endif
00139             enddo ! k = 1, number_of_fields_allocated
00140          enddo ! j = 1, num_blocks
00141       enddo ! i = 1, size (point_data)
00142 
00143    else
00144 
00145       nullify (method_final)
00146    endif ! associated (point_data)
00147 
00148 !===============================================================================
00149 !  2nd merge mask
00150 !===============================================================================
00151 
00152    method_patches(:)%new_mask_id = prism_undefined
00153 
00154    if (associated (mask_data)) then
00155 
00156       allocate (mask_final(size (mask_data)))
00157 
00158       do i = 1, size (mask_data)
00159 
00160          num_blocks = get_num_values (mask_id_map, i)
00161 
00162          ! get the mask_ids of the mask associated to the current user mask id
00163          mask_ids(1:num_blocks) = get_values (mask_id_map, i, num_blocks)
00164 
00165          do j = 1, num_blocks
00166 
00167             k = get_value (mask_id_map, i, j)
00168 
00169             if (j == 1) then
00170 
00171                ! make a copy of the first mask
00172                call copy_mask (mask_final (i), Masks (k))
00173 
00174                ! adjust mask shape in order be able to store mask information of all blocks
00175                call adjust_mask_shape (mask_final (i), mask_ids (1:num_blocks))
00176 
00177             else
00178 
00179                ! merge current field with first field (copy ta)
00180                call merge_mask (mask_final(i), Masks (k))
00181             endif
00182 
00183             ! for all fields
00184             do k = 1, number_of_fields_allocated
00185 
00186                ! find fields that reference the current mask
00187                if (Fields(k)%status /= psmile_status_free .and. &
00188                    Fields(k)%mask_id == mask_ids(j)) then
00189 
00190                   method_patches(k)%new_mask_id = i
00191 
00192                endif
00193             enddo ! k = 1, number_of_fields_allocated
00194          enddo ! j = 1, num_blocks
00195       enddo ! i = 1, size (mask_data)
00196 
00197    else
00198 
00199       nullify (mask_final)
00200    endif ! (associated (mask_data))
00201 
00202 !===============================================================================
00203 !  3rd merge vars/fields
00204 !===============================================================================
00205 
00206    if (associated (var_data)) then
00207 
00208       allocate (field_final(size (var_data)))
00209 
00210       !for all transients
00211       do i = 1, size (var_data)
00212 
00213          num_blocks = get_num_values (var_id_map, i)
00214 
00215          ! get the var_ids of the vars associated to the current user var id
00216          var_ids(1:num_blocks) = get_values (var_id_map, i, num_blocks)
00217 
00218          do j = 1, num_blocks
00219 
00220             k = get_value (var_id_map, i, j)
00221 
00222             if (j == 1) then
00223 
00224                ! make a copy of the first field and its Taskin/out (we do not need to adjust the
00225                ! send/recv_indices in taskin/out for the first field)
00226                call copy_field_information (field_final (i), Fields (k),     &
00227                                             method_patches(k)%new_method_id, &
00228                                             method_patches(k)%new_mask_id)
00229 
00230             else
00231 
00232                ! merge current field with first field (copy ta)
00233                call merge_fields (field_final(i), Fields (k), method_patches(k)%info_offset)
00234             endif
00235 
00236             call add_pair (old_to_new_var_id, k ,i)
00237 
00238          enddo ! j = 1, num_blocks
00239       enddo ! i = 1, size (var_data)
00240    else
00241       nullify (field_final)
00242    endif
00243 
00244 !===============================================================================
00245 !  4th cleanup
00246 !===============================================================================
00247 
00248    ! free old data structures, might be omitted later if the data needs to be reused. One needs
00249    ! to be careful when deleting data, because the new Fields and Methods are in some parts just
00250    ! shallow copies of the original data structures
00251    call psmile_deallocate_fields (ierror)
00252    call psmile_deallocate_methods (ierror)
00253    call psmile_deallocate_masks (ierror)
00254 
00255    ! set data structures to new ones
00256    Fields  => field_final
00257    Methods => method_final
00258    Masks   => mask_final
00259 
00260 !    Number_of_Masks_allocated   = size (Masks)
00261    if (associated (Methods)) then
00262       Number_of_Methods_allocated = size (Methods)
00263    else
00264       Number_of_Methods_allocated = 0
00265    endif
00266    if (associated (Fields)) then
00267       Number_of_Fields_allocated  = size (Fields)
00268    else
00269       Number_of_Fields_allocated = 0
00270    endif
00271    if (associated (Masks)) then
00272       Number_of_Masks_allocated  = size (Masks)
00273    else
00274       Number_of_Masks_allocated = 0
00275    endif
00276 
00277    ! we need to do some adjustments for fields that use userdefined interpolation
00278 
00279    ! for all fields
00280    do i = 1, Number_of_Fields_allocated
00281       call adjust_assoc_var_ids (Fields(i))
00282    enddo
00283 
00284    call free_multimap (old_to_new_var_id)
00285 
00286 #ifdef VERBOSE
00287    print 9980, trim(ch_id), ierror
00288    call psmile_flushstd
00289 #endif /* VERBOSE */
00290 
00291 9990 format (1x, a, ': psmile_merge_fields: ')
00292 9980 format (1x, a, ': psmile_merge_fields: eof ierror =', i5)
00293 
00294 contains
00295 
00296    ! ===========================================================================
00297 
00298    subroutine copy_method (method_to, method_from)
00299 
00300       type (method), intent (in)  :: method_from
00301       type (method), intent (out) :: method_to
00302 
00303       call debug_routine_start("copy_method")
00304 
00305       ! this is just a shallow copy of the method, which means that
00306       ! all pointers of method_to point to the same objects in memory
00307       ! than method_from
00308       method_to = method_from
00309 
00310       ! we need to avoid having multiple pointers to the same objects,
00311       ! because this causes problems in psmile_deallocate
00312       nullify (method_to%coords_pointer, method_to%subgrid_pointer, &
00313                method_to%vector_pointer, method_to%halo_pointer)
00314 
00315       ! because we want to adjust the send and recv infos we have to do a deep copy of them
00316       if (method_from%n_send_info_direct > 0) then
00317 
00318          allocate (method_to%send_infos_direct(method_from%n_send_info_direct))
00319          method_to%send_infos_direct = &
00320             method_from%send_infos_direct(1:method_from%n_send_info_direct)
00321       endif
00322       if (method_from%n_send_info_coupler > 0) then
00323 
00324          allocate (method_to%send_infos_coupler(method_from%n_send_info_coupler))
00325          method_to%send_infos_coupler = &
00326             method_from%send_infos_coupler(1:method_from%n_send_info_coupler)
00327       endif
00328       if (method_from%n_send_info_appl > 0) then
00329 
00330          allocate (method_to%send_infos_appl(method_from%n_send_info_appl))
00331          method_to%send_infos_appl = &
00332             method_from%send_infos_appl(1:method_from%n_send_info_appl)
00333       endif
00334       if (method_from%n_recv_info_direct > 0) then
00335 
00336          allocate (method_to%recv_infos_direct(method_from%n_recv_info_direct))
00337          method_to%recv_infos_direct = &
00338             method_from%recv_infos_direct(1:method_from%n_recv_info_direct)
00339       endif
00340       if (method_from%n_recv_info_coupler > 0) then
00341 
00342          allocate (method_to%recv_infos_coupler(method_from%n_recv_info_coupler))
00343          method_to%recv_infos_coupler = &
00344             method_from%recv_infos_coupler(1:method_from%n_recv_info_coupler)
00345       endif
00346 
00347       call debug_routine_end("copy_method")
00348 
00349    end subroutine copy_method
00350 
00351    ! ===========================================================================
00352 
00353    subroutine merge_method (method_to, method_from, info_offset)
00354 
00355       type (method), intent (in)             :: method_from
00356       type (method), intent (inout)          :: method_to
00357       type (info_index_offset), intent (out) :: info_offset
00358 
00359       call debug_routine_start("merge_method")
00360 
00361       info_offset%send_direct  = method_to%n_send_info_direct
00362       info_offset%send_coupler = method_to%n_send_info_coupler
00363       info_offset%send_appl    = method_to%n_send_info_appl
00364       info_offset%recv_direct  = method_to%n_recv_info_direct
00365       info_offset%recv_coupler = method_to%n_recv_info_coupler
00366 
00367       if (method_from%n_send_info_direct > 0) then
00368 
00369          method_to%send_infos_direct => psmile_realloc(method_to%send_infos_direct, &
00370                                                        method_to%n_send_info_direct + &
00371                                                        method_from%n_send_info_direct)
00372 
00373          method_to%send_infos_direct(method_to%n_send_info_direct+1:) = &
00374             method_from%send_infos_direct(1:method_from%n_send_info_direct)
00375 
00376          method_to%n_send_info_direct = method_to%n_send_info_direct + &
00377                                         method_from%n_send_info_direct
00378       endif
00379 
00380       if (method_from%n_send_info_coupler > 0) then
00381 
00382          method_to%send_infos_coupler => psmile_realloc(method_to%send_infos_coupler, &
00383                                                        method_to%n_send_info_coupler + &
00384                                                        method_from%n_send_info_coupler)
00385 
00386          method_to%send_infos_coupler(method_to%n_send_info_coupler+1:) = &
00387             method_from%send_infos_coupler(1:method_from%n_send_info_coupler)
00388 
00389          method_to%n_send_info_coupler = method_to%n_send_info_coupler + &
00390                                         method_from%n_send_info_coupler
00391       endif
00392 
00393       if (method_from%n_send_info_appl > 0) then
00394 
00395          method_to%send_infos_appl => psmile_realloc(method_to%send_infos_appl, &
00396                                                        method_to%n_send_info_appl + &
00397                                                        method_from%n_send_info_appl)
00398 
00399          method_to%send_infos_appl(method_to%n_send_info_appl+1:) = &
00400             method_from%send_infos_appl(1:method_from%n_send_info_appl)
00401 
00402          method_to%n_send_info_appl = method_to%n_send_info_appl + &
00403                                         method_from%n_send_info_appl
00404       endif
00405 
00406       if (method_from%n_recv_info_direct > 0) then
00407 
00408          method_to%recv_infos_direct => psmile_realloc(method_to%recv_infos_direct, &
00409                                                        method_to%n_recv_info_direct + &
00410                                                        method_from%n_recv_info_direct)
00411 
00412          method_to%recv_infos_direct(method_to%n_recv_info_direct+1:) = &
00413             method_from%recv_infos_direct(1:method_from%n_recv_info_direct)
00414 
00415          method_to%n_recv_info_direct = method_to%n_recv_info_direct + &
00416                                         method_from%n_recv_info_direct
00417       endif
00418 
00419       if (method_from%n_recv_info_coupler > 0) then
00420 
00421          method_to%recv_infos_coupler => psmile_realloc(method_to%recv_infos_coupler, &
00422                                                        method_to%n_recv_info_coupler + &
00423                                                        method_from%n_recv_info_coupler)
00424 
00425          method_to%recv_infos_coupler(method_to%n_recv_info_coupler+1:) = &
00426             method_from%recv_infos_coupler(1:method_from%n_recv_info_coupler)
00427          method_to%n_recv_info_coupler = method_to%n_recv_info_coupler + &
00428                                         method_from%n_recv_info_coupler
00429       endif
00430 
00431       call debug_routine_end("merge_method")
00432 
00433    end subroutine merge_method
00434 
00435    ! ===========================================================================
00436 
00437    subroutine adjust_send_field_information (send_field_infos, offset)
00438 
00439       type (send_field_information), intent (inout) :: send_field_infos(:)
00440       integer, intent (in)                          :: offset
00441 
00442       integer :: i
00443 
00444       do i = 1, size (send_field_infos)
00445          send_field_infos(i)%send_info_index = send_field_infos(i)%send_info_index + offset
00446       enddo
00447 
00448    end subroutine
00449 
00450    ! ===========================================================================
00451 
00452    subroutine adjust_recv_field_information (recv_field_infos, offset)
00453 
00454       type (recv_field_information), intent (inout) :: recv_field_infos(:)
00455       integer, intent (in)                          :: offset
00456 
00457       integer :: i
00458 
00459       do i = 1, size (recv_field_infos)
00460          recv_field_infos(i)%recv_info_index = recv_field_infos(i)%recv_info_index + offset
00461       enddo
00462 
00463    end subroutine
00464 
00465    ! ===========================================================================
00466 
00467    function get_max_num_blocks ()
00468 
00469       integer :: get_max_num_blocks
00470       integer :: i
00471 
00472       get_max_num_blocks = 1
00473 
00474       if (associated (grid_data)) then
00475 
00476          do i = 1, size (grid_data)
00477             get_max_num_blocks = max (get_max_num_blocks, &
00478                                     get_num_values(grid_id_map, i))
00479          enddo
00480       else
00481          get_max_num_blocks = 0
00482       endif
00483 
00484    end function get_max_num_blocks
00485 
00486    ! ===========================================================================
00487 
00488    subroutine insert_sub_array (tgt, tgt_shape, src, src_shape)
00489 
00490       integer, intent (in) :: tgt_shape(2,ndim_3d)
00491 
00492       logical, intent (out)   :: tgt(tgt_shape(1,1):tgt_shape(2,1), 
00493                                      tgt_shape(1,2):tgt_shape(2,2), 
00494                                      tgt_shape(1,3):tgt_shape(2,3))
00495 
00496       integer, intent (in) :: src_shape(2,ndim_3d)
00497 
00498       logical, intent (in)    :: src(src_shape(1,1):src_shape(2,1), 
00499                                      src_shape(1,2):src_shape(2,2), 
00500                                      src_shape(1,3):src_shape(2,3))
00501 
00502       tgt(src_shape(1,1):src_shape(2,1), &
00503           src_shape(1,2):src_shape(2,2), &
00504           src_shape(1,3):src_shape(2,3)) = src
00505 
00506    end subroutine insert_sub_array
00507 
00508    ! ===========================================================================
00509 
00510    subroutine adjust_mask_shape (mask_data, mask_ids)
00511 
00512       type (mask), intent (inout) :: mask_data
00513       integer, intent (in)        :: mask_ids (:)
00514 
00515       integer          :: mask_actual_shape(2, ndim_3d)
00516       integer          :: i, size_of_shape(2)
00517       logical, pointer :: new_mask_array(:)
00518 
00519       size_of_shape = get_size_of_shape (Grids(mask_data%grid_id)%grid_type)
00520 
00521       do i = 1, size_of_shape(2)
00522 
00523          mask_actual_shape(1, i) = minval (Masks(mask_ids)%mask_shape(1, i))
00524          mask_actual_shape(2, i) = maxval (Masks(mask_ids)%mask_shape(2, i))
00525       enddo
00526 
00527       do i = size_of_shape(2) + 1, ndim_3d
00528 
00529          mask_actual_shape(1, i) = 1
00530          mask_actual_shape(2, i) = 1
00531       enddo
00532 
00533       nullify (new_mask_array)
00534       new_mask_array => psmile_realloc(new_mask_array, product (mask_actual_shape(2,:) - &
00535                                                                 mask_actual_shape(1,:) + 1))
00536 
00537       call insert_sub_array (new_mask_array, mask_actual_shape, &
00538                              mask_data%mask_array, mask_data%mask_shape)
00539 
00540       mask_data%mask_shape = mask_actual_shape
00541       mask_data%mask_array => new_mask_array
00542 
00543    end subroutine adjust_mask_shape
00544 
00545    ! ===========================================================================
00546 
00547    subroutine copy_mask (mask_to, mask_from)
00548 
00549       type (mask), intent (in)  :: mask_from
00550       type (mask), intent (out) :: mask_to
00551 
00552       logical, pointer :: mask_array
00553 
00554       call debug_routine_start("copy_mask")
00555 
00556       mask_to = mask_from
00557 
00558       nullify (mask_to%mask_array)
00559       mask_to%mask_array => psmile_realloc (mask_to%mask_array, size (mask_from%mask_array))
00560 
00561       mask_to%mask_array = mask_from%mask_array
00562 
00563       call debug_routine_end("copy_mask")
00564 
00565    end subroutine copy_mask
00566 
00567    ! ===========================================================================
00568 
00569    subroutine merge_mask (mask_to, mask_from)
00570 
00571       type (mask), intent (in)    :: mask_from
00572       type (mask), intent (inout) :: mask_to
00573 
00574       call debug_routine_start("merge_mask")
00575 
00576       call insert_sub_array (mask_to%mask_array, mask_to%mask_shape, &
00577                              mask_from%mask_array, mask_from%mask_shape)
00578 
00579       call debug_routine_end("merge_mask")
00580 
00581    end subroutine merge_mask
00582 
00583    ! ===========================================================================
00584 
00585    subroutine copy_taskin_type (taskin_to, taskin_from)
00586 
00587       type (taskin_type), intent (in)  :: taskin_from
00588       type (taskin_type), intent (out) :: taskin_to
00589 
00590       call debug_routine_start("copy_taskin_type")
00591 
00592       ! this is just a shallow copy of the taskin, which means that
00593       ! all pointers of taskin_to point to the same objects in memory
00594       ! than taskin_from
00595       taskin_to = taskin_from
00596 
00597       ! make a deep copy of taskin_from%In_channel
00598       if (associated (taskin_from%In_channel)) then
00599          nullify (taskin_to%In_channel)
00600          taskin_to%In_channel => psmile_realloc (taskin_to%In_channel, &
00601                                                  size (taskin_from%In_channel))
00602          taskin_to%In_channel = taskin_from%In_channel
00603       endif
00604 
00605       if (taskin_from%n_recv_direct > 0) then
00606 
00607          allocate (taskin_to%recv_direct(taskin_from%n_recv_direct))
00608          taskin_to%recv_direct = taskin_from%recv_direct(1:taskin_from%n_recv_direct)
00609          taskin_to%n_recv_direct = taskin_from%n_recv_direct
00610       else
00611          taskin_to%n_recv_direct = 0
00612          nullify (taskin_to%recv_direct)
00613       endif
00614 
00615       if (taskin_from%n_recv_coupler > 0) then
00616 
00617          allocate (taskin_to%recv_coupler(taskin_from%n_recv_coupler))
00618          taskin_to%recv_coupler = taskin_from%recv_coupler(1:taskin_from%n_recv_coupler)
00619          taskin_to%n_recv_coupler = taskin_from%n_recv_coupler
00620       else
00621          nullify (taskin_to%recv_coupler)
00622          taskin_to%n_recv_coupler = 0
00623       endif
00624 
00625       call debug_routine_end("copy_taskin_type")
00626 
00627    end subroutine copy_taskin_type
00628 
00629    ! ===========================================================================
00630 
00631    subroutine copy_taskout_type (taskout_to, taskout_from)
00632 
00633       type (taskout_type), intent (in)  :: taskout_from(:)
00634       type (taskout_type), intent (out) :: taskout_to(:)
00635 
00636       integer :: i
00637 
00638       call debug_routine_start("copy_taskout_type")
00639 
00640       ! this is just a shallow copy of the taskout, which means that
00641       ! all pointers of taskout_to point to the same objects in memory
00642       ! than taskout_from
00643       taskout_to = taskout_from
00644 
00645       do i = 1, size (taskout_from)
00646 
00647          if (taskout_from(i)%n_send_direct > 0) then
00648 
00649             allocate (taskout_to(i)%send_direct(taskout_from(i)%n_send_direct))
00650             taskout_to(i)%send_direct = taskout_from(i)%send_direct(1:taskout_from(i)%n_send_direct)
00651             taskout_to(i)%n_send_direct = taskout_from(i)%n_send_direct
00652          else
00653             nullify (taskout_to(i)%send_direct)
00654             taskout_to(i)%n_send_direct = 0
00655          endif
00656 
00657          if (taskout_from(i)%n_send_coupler > 0) then
00658 
00659             allocate (taskout_to(i)%send_coupler(taskout_from(i)%n_send_coupler))
00660             taskout_to(i)%send_coupler = taskout_from(i)%send_coupler(1:taskout_from(i)%n_send_coupler)
00661             taskout_to(i)%n_send_coupler = taskout_from(i)%n_send_coupler
00662          else
00663             nullify (taskout_to(i)%send_coupler)
00664             taskout_to(i)%n_send_coupler = 0
00665          endif
00666 
00667          if (taskout_from(i)%n_send_appl > 0) then
00668 
00669             allocate (taskout_to(i)%send_appl(taskout_from(i)%n_send_appl))
00670             taskout_to(i)%send_appl = taskout_from(i)%send_appl(1:taskout_from(i)%n_send_appl)
00671             taskout_to(i)%n_send_appl = taskout_from(i)%n_send_appl
00672          else
00673             nullify (taskout_to(i)%send_appl)
00674             taskout_to(i)%n_send_appl = 0
00675          endif
00676       enddo
00677 
00678       call debug_routine_end("copy_taskout_type")
00679 
00680    end subroutine copy_taskout_type
00681 
00682    ! ===========================================================================
00683 
00684    subroutine copy_field_information (field_to, field_from, new_method_id, new_mask_id)
00685 
00686       type (gridfunction), intent (in)  :: field_from
00687       type (gridfunction), intent (out) :: field_to
00688       integer, intent (in)              :: new_method_id
00689       integer, intent (in)              :: new_mask_id
00690 
00691       call debug_routine_start("copy_field_information")
00692 
00693       field_to = field_from
00694 
00695       field_to%method_id = new_method_id
00696       field_to%mask_id   = new_mask_id
00697 
00698       call copy_taskin_type(field_to%Taskin, field_from%Taskin)
00699 
00700       if (associated (field_from%Taskout)) then
00701          if (size (field_from%Taskout) > 0) then
00702 
00703             allocate (field_to%Taskout(size (field_from%Taskout)))
00704             call copy_taskout_type(field_to%Taskout, field_from%Taskout)
00705 
00706          endif
00707       endif
00708 
00709       ! I do not know what to do with the io stuff, therefore I am leaving this for later
00710 
00711       call debug_routine_end("copy_field_information")
00712 
00713    end subroutine copy_field_information
00714 
00715    ! ===========================================================================
00716 
00717    subroutine merge_taskin (taskin_to, taskin_from, info_offset)
00718 
00719       type (taskin_type), intent (in)       :: taskin_from
00720       type (taskin_type), intent (inout)    :: taskin_to
00721       type (info_index_offset), intent (in) :: info_offset
00722 
00723       call debug_routine_start("merge_taskin")
00724 
00725       if (taskin_from%n_recv_direct > 0) then
00726 
00727          taskin_to%recv_direct => psmile_realloc (taskin_to%recv_direct, &
00728                                                   taskin_to%n_recv_direct + &
00729                                                   taskin_from%n_recv_direct)
00730 
00731          taskin_to%recv_direct(taskin_to%n_recv_direct+1:) = &
00732             taskin_from%recv_direct(1:taskin_from%n_recv_direct)
00733 
00734          call adjust_recv_field_information (taskin_to%recv_direct(taskin_to%n_recv_direct+1:), &
00735                                              info_offset%recv_direct)
00736          taskin_to%n_recv_direct = taskin_to%n_recv_direct + taskin_from%n_recv_direct
00737       endif
00738 
00739       if (taskin_from%n_recv_coupler > 0) then
00740 
00741          taskin_to%recv_coupler => psmile_realloc (taskin_to%recv_coupler, &
00742                                                   taskin_to%n_recv_coupler + &
00743                                                   taskin_from%n_recv_coupler)
00744 
00745          taskin_to%recv_coupler(taskin_to%n_recv_coupler+1:) = &
00746             taskin_from%recv_coupler(1:taskin_from%n_recv_coupler)
00747 
00748          call adjust_recv_field_information (taskin_to%recv_coupler(taskin_to%n_recv_coupler+1:), &
00749                                              info_offset%recv_coupler)
00750 
00751          taskin_to%n_recv_coupler = taskin_to%n_recv_coupler + taskin_from%n_recv_coupler
00752       endif
00753 
00754       call debug_routine_end("merge_taskin")
00755 
00756    end subroutine merge_taskin
00757 
00758    ! ===========================================================================
00759 
00760    subroutine merge_taskout (taskout_to, taskout_from, info_offset)
00761 
00762       type (taskout_type), intent (in)      :: taskout_from(:)
00763       type (taskout_type), intent (inout)   :: taskout_to(:)
00764       type (info_index_offset), intent (in) :: info_offset
00765 
00766       integer :: i
00767 
00768       call debug_routine_start("merge_taskout")
00769 
00770 #ifdef PRISM_ASSERTION
00771       if (size (taskout_to) /= size (taskout_from)) then
00772          call psmile_assert (__FILE__, __LINE__, &
00773                              "merge_taskout taskout_to and taskout_from have different sizes")
00774       endif
00775 #endif
00776 
00777       do i = 1, size (taskout_to)
00778 
00779          if (taskout_from(i)%n_send_direct > 0) then
00780 
00781             taskout_to(i)%send_direct => psmile_realloc (taskout_to(i)%send_direct, &
00782                                                          taskout_to(i)%n_send_direct + &
00783                                                          taskout_from(i)%n_send_direct)
00784 
00785             taskout_to(i)%send_direct(taskout_to(i)%n_send_direct+1:) = &
00786                taskout_from(i)%send_direct(1:taskout_from(i)%n_send_direct)
00787 
00788             call adjust_send_field_information (taskout_to(i)%send_direct(taskout_to(i)%n_send_direct+1:), &
00789                                                 info_offset%send_direct)
00790             taskout_to(i)%n_send_direct = taskout_to(i)%n_send_direct + taskout_from(i)%n_send_direct
00791          endif
00792 
00793          if (taskout_from(i)%n_send_coupler > 0) then
00794 
00795             taskout_to(i)%send_coupler => psmile_realloc (taskout_to(i)%send_coupler, &
00796                                                           taskout_to(i)%n_send_coupler + &
00797                                                           taskout_from(i)%n_send_coupler)
00798 
00799             taskout_to(i)%send_coupler(taskout_to(i)%n_send_coupler+1:) = &
00800                taskout_from(i)%send_coupler(1:taskout_from(i)%n_send_coupler)
00801 
00802             call adjust_send_field_information (taskout_to(i)%send_coupler(taskout_to(i)%n_send_coupler+1:), &
00803                                                 info_offset%send_coupler)
00804             taskout_to(i)%n_send_coupler = taskout_to(i)%n_send_coupler + taskout_from(i)%n_send_coupler
00805          endif
00806 
00807          if (taskout_from(i)%n_send_appl > 0) then
00808 
00809             taskout_to(i)%send_appl => psmile_realloc (taskout_to(i)%send_appl, &
00810                                                        taskout_to(i)%n_send_appl + &
00811                                                        taskout_from(i)%n_send_appl)
00812 
00813             taskout_to(i)%send_appl(taskout_to(i)%n_send_appl+1:) = &
00814                taskout_from(i)%send_appl(1:taskout_from(i)%n_send_appl)
00815 
00816             call adjust_send_field_information (taskout_to(i)%send_appl(taskout_to(i)%n_send_appl+1:), &
00817                                                 info_offset%send_appl)
00818             taskout_to(i)%n_send_appl = taskout_to(i)%n_send_appl + taskout_from(i)%n_send_appl
00819          endif
00820       enddo ! i = 1, size (taskout_to)
00821 
00822       call debug_routine_end("merge_taskout")
00823 
00824    end subroutine merge_taskout
00825 
00826    ! ===========================================================================
00827 
00828    subroutine merge_fields (field_to, field_from, info_offset)
00829 
00830       type (gridfunction), intent (in)      :: field_from
00831       type (gridfunction), intent (inout)   :: field_to
00832       type (info_index_offset), intent (in) :: info_offset
00833 
00834       call debug_routine_start("merge_fields")
00835 
00836       ! I do not know what to do with IO stuff...but it did not work for multiblock before as well...
00837 
00838       call merge_taskin (field_to%taskin, field_from%taskin, info_offset)
00839       if (associated (field_from%taskout)) then
00840          call merge_taskout (field_to%taskout, field_from%taskout, info_offset)
00841       else
00842          nullify (field_to%taskout)
00843       endif
00844 
00845       call debug_routine_end("merge_fields")
00846 
00847    end subroutine merge_fields
00848 
00849    ! ===========================================================================
00850 
00851    subroutine adjust_assoc_var_ids (field)
00852 
00853       type (gridfunction), intent (inout) :: field
00854 
00855       integer :: i
00856 
00857       call debug_routine_start ("adjust_assoc_var_ids")
00858 
00859       ! if there taskout's defined for this field
00860       if (associated (field%taskout)) then
00861 
00862          ! for all taskout's
00863          do i = 1, size (field%taskout)
00864 
00865             ! if the current taskout has a associated var_id
00866             if (field%taskout(i)%assoc_var_id /= psmile_undef) then
00867 
00868                ! if we have a new var id for the current assoc_var_id
00869                if (is_valid_key (old_to_new_var_id, field%taskout(i)%assoc_var_id)) then
00870 
00871                   ! adjust assoc_var_id
00872                   field%taskout(i)%assoc_var_id = get_value (old_to_new_var_id, &
00873                                                              field%taskout(i)%assoc_var_id, 1)
00874                endif
00875             endif
00876          enddo ! i
00877       endif ! associated (field%taskout)
00878 
00879       ! for all inchannels of the field
00880       do i = 1, field%taskin%nbr_inchannels
00881 
00882          ! if the current inchannel has a associated var_id
00883          if (field%taskin%in_channel(i)%assoc_var_id /= psmile_undef) then
00884 
00885             ! if we have a new var id for the current assoc_var_id
00886             if (is_valid_key (old_to_new_var_id, field%taskin%in_channel(i)%assoc_var_id)) then
00887 
00888                ! adjust assoc_var_id
00889                field%taskin%in_channel(i)%assoc_var_id = get_value (old_to_new_var_id, &
00890                                                       field%taskin%in_channel(i)%assoc_var_id, 1)
00891             endif
00892          endif
00893       enddo ! i
00894 
00895       call debug_routine_end ("adjust_assoc_var_ids")
00896 
00897    end subroutine adjust_assoc_var_ids
00898 
00899    ! ===========================================================================
00900 
00901    subroutine debug_routine_start (name)
00902 
00903       character (len=*), intent (in) :: name
00904 #ifdef VERBOSE
00905       print "(1x, a, ': ', a,': ')", trim(ch_id), name
00906       call psmile_flushstd
00907 #endif /* VERBOSE */
00908    end subroutine debug_routine_start
00909 
00910    ! ===========================================================================
00911 
00912    subroutine debug_routine_end (name)
00913 
00914       character (len=*), intent (in) :: name
00915 #ifdef VERBOSE
00916       print "(1x, a, ': ', a,': eof')", trim(ch_id), name
00917       call psmile_flushstd
00918 #endif /* VERBOSE */
00919    end subroutine debug_routine_end
00920 
00921 end subroutine psmile_merge_fields

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1