boundary.f90

Path: src/util/boundary.f90
Last Update: Sat Apr 23 00:01:52 JST 2005

    Copyright (C) GFD Dennou Club, 2004. All rights reserved.

begin

Subroutine Boundary

  * Developer: SUGIYAMA Ko-ichiro (sugiyama@gfd-dennou.org)
  * Version: $Id: boundary.f90,v 1.2 2005/04/22 15:01:52 sugiyama Exp $
  * Tag Name: $Name:  $
  * Change History:

Overview

境界条件を与えるためのモジュール. 引数 type に指定された境界条件から, 「のり代」の部分の値を設定する

Error Handling

Known Bugs

Note

水平鉛直に摩擦無し境界とする設定は, 移流のテスト計算用であることに注意.

Future Plans

end

Methods

boundary  

Included Modules

dc_trace dc_message gridset

Public Instance methods

type :character(*), intent(in)
: end — 暗黙の型宣言禁止 begin
 Input
var(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(inout)
: In/Out

[Source]

subroutine boundary(type, var)
                                                                 !=begin
  !== Dependency

                                                                 !=end  
  !--- 暗黙の型宣言禁止
  implicit none
                                                                 !=begin  
  !== Input
  character(*), intent(in)   :: type

  !== In/Out
  real(8), intent(inout)     :: var(DimXMin:DimXMax, DimZMin:DimZMax)
                                                                 !=end    
  !--- 内部変数
  integer   :: i

   call BeginSub("boundary",                  fmt="%c",                        c1="Set boundary conditions.")

  if (type == "CC") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do
     
     !--- z 方向: 周期境界     
     do i = 0, Margin
        var(:, RegZMin - i) = var(:, RegZMax - i) 
     end do
     do i = 1, Margin
        var(:, RegZMax + i) = var(:, RegZMin + i) 
     end do

  elseif (type == "CsS") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do

     !--- z 方向: 対称
     do i = 0, Margin
        var(:, RegZMin - i) = var(:, RegZMin + 1 + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax + 1 - i )
     end do
     
  elseif (type == "CfS") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do

     !--- z 方向: 対称
     var(:, RegZMin) = 0.0d0
     var(:, RegZMax) = 0.0d0
     do i = 1, Margin
        var(:, RegZMin - i) = var(:, RegZMin + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax - i )
     end do
     
  elseif (type == "CsA") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do

     !--- z 方向: 反対称     
     do i = 0, Margin
        var(:, RegZMin - i) = - var(:, RegZMin + 1 + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i) = - var(:, RegZMax + 1 - i )
     end do
     
  elseif (type == "CfA") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do
     
     !--- z 方向: 反対称     
     var(:, RegZMin) = 0.0d0
     var(:, RegZMax) = 0.0d0
     do i = 1, Margin
        var(:, RegZMin - i) = - var(:, RegZMin + i )        
    end do
     do i = 1, Margin
        var(:, RegZMax + i) = - var(:, RegZMax - i )
     end do

     
  !--- test     
  elseif (type == "CsAS") then 
     !--- x 方向: 周期境界
     do i = 0, Margin
        var(RegXMin - i , RegZMin:RegZMax)  = var(RegXMax - i , RegZMin:RegZMax) 
     end do
     do i = 1, Margin
        var(RegXMax + i , RegZMin:RegZMax)  = var(RegXMin + i , RegZMin:RegZMax) 
     end do

     !--- z 方向: 対称と反対称     
     do i = 0, Margin
        var(:, RegZMin - i) = - var(:, RegZMin + 1 + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax + 1 - i )
     end do


  elseif (type == "sSsS") then 
     !--- x 方向: 対称
     do i = 0, Margin
        var(RegXMin - i,:) = var(RegXMin + 1 + i,: )
     end do
     do i = 1, Margin
        var(RegXMax + i ,:) = var(RegXMax + 1 - i,: )
     end do

     !--- z 方向: 対称
     do i = 0, Margin
        var(:, RegZMin - i) = var(:, RegZMin + 1 + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax + 1 - i )
     end do
     

  elseif (type == "sSfS") then 
     !--- x 方向: 対称
     do i = 0, Margin
        var(RegXMin - i,:) = var(RegXMin + 1 + i,: )
     end do
     do i = 1, Margin
        var(RegXMax + i ,:) = var(RegXMax + 1 - i,: )
     end do

     !--- z 方向: 対称
     var(:,RegZMin) = var(:,RegZMin + 1)

     do i = 1, Margin
        var(:, RegZMin - i) = var(:, RegZMin + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax - i )
     end do

  elseif (type == "fSsS") then 
     !--- x 方向: 対称
     var(RegXMin,:) = var(RegXMin + 1,:)

     do i = 1, Margin
        var(RegXMin - i,:) = var(RegXMin + i,: )
     end do
     do i = 1, Margin
        var(RegXMax + i ,:) = var(RegXMax - i,: )
     end do

     !--- z 方向: 対称
     do i = 0, Margin
        var(:, RegZMin - i) = var(:, RegZMin + 1 + i )
     end do
     do i = 1, Margin
        var(:, RegZMax + i ) = var(:, RegZMax + 1 - i )
     end do
     

  else
     call MessageNotify("Error", "boundary",                        "Unknown Boundary Type", c1=type)

  end if

  call EndSub("boundary")

end subroutine boundary

[Validate]