historyautoputaxismpi.f90

Path: gtool/gtool_historyauto/historyautoputaxismpi.f90
Last Update: Sun May 10 21:19:17 +0900 2009

領域全体の座標データの設定 (MPI 用)

Settings of data on whole area of axes (for MPI)

Authors:Yasuhiro MORIKAWA
Version:$Id: historyautoputaxismpi.f90,v 1.1 2009-05-10 12:19:17 morikawa Exp $
Tag Name:$Name: gtool5-20101228-1 $
Copyright:Copyright (C) GFD Dennou Club, 2008-2009. All rights reserved.
License:See COPYRIGHT

Required files

Methods

Included Modules

gtool_historyauto_internal gtool_history dc_trace dc_error dc_message dc_types

Public Instance methods

Subroutine :
dim :character(*), intent(in)
: 座標変数の名称.

ここで指定するものは, HistoryAutoCreate の 引数 dims で既に指定されてい なければなりません.

Name of dimensional variable.

This name must be specified by an argument dims in "HistoryAutoCreate".

array(:) :real(DP), intent(in)
: 座標データ.

Data of axes.

MPI 使用時に, 各々のノード上のデータを単一ファイルに 集約して出力する場合には, このサブルーチンに領域全体の座標データを与えてください. また, HistoryAutoCreate のオプショナル論理型引数 flag_mpi_gather に .true. を与えてください.

When MPI is used, if data on each node is integrated and output to one file, give data of axes in whole area to this subroutine. And give .true. to optional logical argument flag_mpi_gather in "HistoryAutoCreate".

[Source]

  subroutine HistoryAutoPutAxisMPIDouble( dim, array )
    !
    ! MPI 使用時に, 各々のノード上のデータを単一ファイルに
    ! 集約して出力する場合には, 
    ! このサブルーチンに領域全体の座標データを与えてください. 
    ! また, HistoryAutoCreate のオプショナル論理型引数 *flag_mpi_gather* 
    ! に .true. を与えてください. 
    !
    ! When MPI is used, if data on each node is integrated and 
    ! output to one file, give data of axes in whole area to 
    ! this subroutine. 
    ! And give .true. to optional logical argument *flag_mpi_gather* 
    ! in "HistoryAutoCreate". 
    !
    use gtool_historyauto_internal, only: initialized, numdims, gthst_axes, data_axes_whole
    use gtool_history, only: HistoryAxisInquire
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, HST_ENOAXISNAME
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, GT_EBADDIMNAME
    use dc_message, only: MessageNotify
    use dc_types, only: DP, STRING, TOKEN
    implicit none
    character(*), intent(in):: dim
                                                  ! 座標変数の名称. 
                              !
                              ! ここで指定するものは, HistoryAutoCreate の
                              ! 引数 *dims* で既に指定されてい
                              ! なければなりません.
                              !
                              ! Name of dimensional variable. 
                              !
                              ! This name must be specified by 
                              ! an argument *dims* in "HistoryAutoCreate". 
                              !
                    
    real(DP), intent(in):: array(:)
                                                  ! 座標データ. 
                              !
                              ! Data of axes. 
                    
    integer:: i, dimsize
    character(STRING):: name
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoPutAxisMPIDouble"
  continue
    call BeginSub(subname, 'dim=%c', c1 = trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name )            ! (out)
      if ( trim(dim) == trim(name) ) then
        dimsize = size( array )
        allocate( data_axes_whole(i) % a_axis( dimsize ) )
        data_axes_whole(i) % a_axis(:) = array(:)
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim

    ! 終了処理, 例外処理
    ! Termination and Exception handling
    !
999 continue
    call StoreError( stat, subname, cause_c = cause_c )
    call EndSub(subname)
  end subroutine HistoryAutoPutAxisMPIDouble
Subroutine :
dim :character(*), intent(in)
array(:) :integer, intent(in)

MPI 使用時に, 各々のノード上のデータを単一ファイルに 集約して出力する場合には, このサブルーチンに領域全体の座標データを与えてください. また, HistoryAutoCreate のオプショナル論理型引数 flag_mpi_gather に .true. を与えてください.

When MPI is used, if data on each node is integrated and output to one file, give data of axes in whole area to this subroutine. And give .true. to optional logical argument flag_mpi_gather in "HistoryAutoCreate".

[Source]

  subroutine HistoryAutoPutAxisMPIInt( dim, array )
    !
    ! MPI 使用時に, 各々のノード上のデータを単一ファイルに
    ! 集約して出力する場合には, 
    ! このサブルーチンに領域全体の座標データを与えてください. 
    ! また, HistoryAutoCreate のオプショナル論理型引数 *flag_mpi_gather* 
    ! に .true. を与えてください. 
    !
    ! When MPI is used, if data on each node is integrated and 
    ! output to one file, give data of axes in whole area to 
    ! this subroutine. 
    ! And give .true. to optional logical argument *flag_mpi_gather* 
    ! in "HistoryAutoCreate". 
    !
    use gtool_historyauto_internal, only: initialized, numdims, gthst_axes, data_axes_whole
    use gtool_history, only: HistoryAxisInquire
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, HST_ENOAXISNAME
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, GT_EBADDIMNAME
    use dc_message, only: MessageNotify
    use dc_types, only: DP, STRING, TOKEN
    implicit none
    character(*), intent(in):: dim
                    
    integer, intent(in):: array(:)
                    
    integer:: i, dimsize
    character(STRING):: name
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoPutAxisMPIInt"
  continue
    call BeginSub(subname, 'dim=%c', c1 = trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name )            ! (out)
      if ( trim(dim) == trim(name) ) then
        dimsize = size( array )
        allocate( data_axes_whole(i) % a_axis( dimsize ) )
        data_axes_whole(i) % a_axis(:) = array(:)
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim

    ! 終了処理, 例外処理
    ! Termination and Exception handling
    !
999 continue
    call StoreError( stat, subname, cause_c = cause_c )
    call EndSub(subname)
  end subroutine HistoryAutoPutAxisMPIInt
Subroutine :
dim :character(*), intent(in)
array(:) :real, intent(in)

MPI 使用時に, 各々のノード上のデータを単一ファイルに 集約して出力する場合には, このサブルーチンに領域全体の座標データを与えてください. また, HistoryAutoCreate のオプショナル論理型引数 flag_mpi_gather に .true. を与えてください.

When MPI is used, if data on each node is integrated and output to one file, give data of axes in whole area to this subroutine. And give .true. to optional logical argument flag_mpi_gather in "HistoryAutoCreate".

[Source]

  subroutine HistoryAutoPutAxisMPIReal( dim, array )
    !
    ! MPI 使用時に, 各々のノード上のデータを単一ファイルに
    ! 集約して出力する場合には, 
    ! このサブルーチンに領域全体の座標データを与えてください. 
    ! また, HistoryAutoCreate のオプショナル論理型引数 *flag_mpi_gather* 
    ! に .true. を与えてください. 
    !
    ! When MPI is used, if data on each node is integrated and 
    ! output to one file, give data of axes in whole area to 
    ! this subroutine. 
    ! And give .true. to optional logical argument *flag_mpi_gather* 
    ! in "HistoryAutoCreate". 
    !
    use gtool_historyauto_internal, only: initialized, numdims, gthst_axes, data_axes_whole
    use gtool_history, only: HistoryAxisInquire
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, HST_ENOAXISNAME
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, GT_EBADDIMNAME
    use dc_message, only: MessageNotify
    use dc_types, only: DP, STRING, TOKEN
    implicit none
    character(*), intent(in):: dim
                    
    real, intent(in):: array(:)
                    
    integer:: i, dimsize
    character(STRING):: name
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoPutAxisMPIReal"
  continue
    call BeginSub(subname, 'dim=%c', c1 = trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name )            ! (out)
      if ( trim(dim) == trim(name) ) then
        dimsize = size( array )
        allocate( data_axes_whole(i) % a_axis( dimsize ) )
        data_axes_whole(i) % a_axis(:) = array(:)
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim

    ! 終了処理, 例外処理
    ! Termination and Exception handling
    !
999 continue
    call StoreError( stat, subname, cause_c = cause_c )
    call EndSub(subname)
  end subroutine HistoryAutoPutAxisMPIReal