!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2004. All rights reserved.
!---------------------------------------------------------------------
                                                                 !=begin
!= Program ArareInit
!
!   * Developer: SUGIYAMA Ko-ichiro
!   * Version: $Id: init.f90,v 1.2.2.8 2006/01/01 08:26:36 kitamo Exp $ 
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview 
!
!$BHs@ENO3X%b%G%k(B deepconv/arare $B$N=i4|CM$r:n@.$9$k(B. 
!$B=i4|CM%U%!%$%k$K$O(B 2 $B%9%F%C%WJ,$N%G!<%?$,I,MW$J$N$G(B, 
!$B:G=i$N(B 1 $B%9%F%C%W$r%*%$%i!<K!$G7W;;$9$k(B
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
! *$BJ}Dx<07O$O=`05=L7O(B.
!
!== Future Plans
!
                                                                 !=end

program ArareInit
                                                                 !=begin
  !== Dependency
  use dc_trace,   only: BeginSub, EndSub, DbgMessage
  use dc_message
  use debugset,   only: debugset_init, DebugOn
  use gridset,    only: gridset_init, DimXMin, DimXMax, DimZMin, DimZMax, &
       &                RegZMin, RegZMax
  use fileset,    only: fileset_init, cfgfile, ReStartFile
  use timeset,    only: timeset_init, NstepShort, DelTimeLong, DelTimeShort
  use nameset,    only: nameset_init
  use physset,    only: physset_init  
  use cloudset,   only: cloudset_init  
  use linlib,     only: linlib_init
  use bcset,      only: bcset_init, ss_BC
  use arareset,   only: arareset_init
  use basicset,   only: basicset_calc_init
  use disturbset, only: disturbset_calc_init, &
       &                    ss_ExnerDisturb, ss_PotTempDisturb, &
       &                    fs_VelXDisturb,  sf_VelZDisturb
  use spongelayer, only: fs_SLayer_fs, sf_SLayer_sf 
                                                                 !=end
  implicit none

  real(8), allocatable  :: ss_Exner_as(:,:)
  real(8), allocatable  :: ss_PotTemp_as(:,:)
  real(8), allocatable  :: ss_DensCloud_as(:,:)
  real(8), allocatable  :: fs_VelX_as(:,:)
  real(8), allocatable  :: sf_VelZ_as(:,:)
  real(8), allocatable  :: ss_Exner_ns(:,:)
  real(8), allocatable  :: ss_PotTemp_ns(:,:)
  real(8), allocatable  :: ss_DensCloud_ns(:,:)
  real(8), allocatable  :: fs_VelX_ns(:,:)
  real(8), allocatable  :: sf_VelZ_ns(:,:)

  real(8), allocatable  :: fs_AdvDiffX_nl(:,:)
  real(8), allocatable  :: sf_AdvDiffZ_nl(:,:)
  real(8), allocatable  :: ss_ExnerTendCond_ns(:,:)
  real(8), allocatable  :: fs_VelXTendAdv_nl(:,:)  ! $B0\N.$K$h$k?eJ?B.EYJQ2=(B
  real(8), allocatable  :: fs_VelXTendTurb_nl(:,:) ! $BMpN.3H;6$K$h$k?eJ?B.EYJQ2=(B
  real(8), allocatable  :: fs_VelXTendNumDiff_nl(:,:) ! $B?tCM3H;6$K$h$k?eJ?B.EYJQ2=(B
  real(8), allocatable  :: sf_VelZTendAdv_nl(:,:)  ! $B0\N.$K$h$k1tD>B.EYJQ2=(B
  real(8), allocatable  :: sf_VelZTendTurb_nl(:,:) ! $BMpN.3H;6$K$h$k1tD>B.EYJQ2=(B
  real(8), allocatable  :: sf_VelZTendBuoy_nl(:,:) ! $BIbNO$K$h$k1tD>B.EYJQ2=(B
  real(8), allocatable  :: sf_VelZTendNumDiff_nl(:,:) ! $B?tCM3H;6$K$h$k1tD>B.EYJQ2=(B  
  real(8), allocatable  :: sf_VelZTendPressGrad_ns(:,:)
                                      ! $B05NO8{G[$K$h$k1tD>B.EYJQ2=(B

  real(8), allocatable  :: fs_VelX_nl(:,:)
  real(8), allocatable  :: sf_VelZ_nl(:,:)
  real(8), allocatable  :: fs_VelX_al(:,:)
  real(8), allocatable  :: sf_VelZ_al(:,:)

  real(8), allocatable  :: ss_PotTemp_nl(:,:)
  real(8), allocatable  :: ss_Exner_nl(:,:)
  real(8), allocatable  :: ss_PotTemp_al(:,:)
  real(8), allocatable  :: ss_Exner_al(:,:)

  real(8), allocatable  :: ss_Km_al(:,:)
  real(8), allocatable  :: ss_Kh_al(:,:)
  real(8), allocatable  :: ss_Km_nl(:,:)
  real(8), allocatable  :: ss_Kh_nl(:,:)

  real(8), allocatable  :: ss_LatHeatPerMass_nl(:,:) ! $BC10L<ANL$"$?$j$N@xG.(B
  real(8), allocatable  :: ss_MassCond_ns(:,:) ! $B6E7k(B/$B>xH/NL(B

  real(8), allocatable  :: ss_DensCloud_nl(:,:) ! $B1@$NL)EY(B(t=0)
  real(8), allocatable  :: ss_DensCloud_al(:,:) ! $B1@$NL)EY(B(t=$B&$(Bt)

  real(8), allocatable  :: ss_SupSatRatio_nl(:,:) ! $B1@$NL)EY(B(t=0)
  real(8), allocatable  :: ss_SupSatRatio_al(:,:) ! $B1@$NL)EY(B(t=$B&$(Bt)
  real(8), allocatable  :: ss_SupSatRatio_ns(:,:) ! $B1@$NL)EY(B(t=0)
  real(8), allocatable  :: ss_SupSatRatio_as(:,:) ! $B1@$NL)EY(B(t=$B&$&S(B)

  real(8), allocatable  :: ss_PotTempTendAdv_nl(:,:)  ! $B0\N.$K$h$k290LJQ2=(B
  real(8), allocatable  :: ss_PotTempTendSrc_nl(:,:)  ! $B@8@.$K$h$k290LJQ2=(B
  real(8), allocatable  :: ss_PotTempTendTurb_nl(:,:) ! $BMpN.3H;6$K$h$k290LJQ2=(B
  real(8), allocatable  :: ss_PotTempTendNumDiff_nl(:,:) ! $B?tCM3H;6$K$h$k290LJQ2=(B

  real(8), allocatable  :: ss_QCond(:,:)        ! $B6E7kG.(B

  real(8)               :: ReStartTime(2)
  integer               :: tau


!!!
!!!$B=i4|2=(B
!!!
  !--- I/O $B%U%!%$%kL>$N=i4|2=(B
  call fileset_init

  !--- $B%G%P%0@_Dj$N=i4|2=(B
  call debugset_init(cfgfile)

  !--- $B3J;RE@>pJs$N=i4|2=(B
  call gridset_init(cfgfile)
  
  !--- $B;~4V@QJ,$N@_Dj$N=i4|2=(B
  call timeset_init(cfgfile)
  
  !--- $B<B83L>$N=i4|2=(B
  call nameset_init(cfgfile)
  
  !--- $BJ*M}NL$N=i4|2=(B
  call physset_init(cfgfile)
  call cloudset_init(cfgfile)
  
  !--- $B7W;;$N@_Dj$rFI$_9~$`(B
  call arareset_init(cfgfile)

  !--- $B6-3&>r7o$r@_Dj(B
  call bcset_init(cfgfile)

  !--- $B<B83@_Dj$rFI$_9~$`(B
  call basicset_calc_init(cfgfile)

  !--- $B<B83@_Dj$rFI$_9~$`(B
  call disturbset_calc_init(cfgfile)
  
  !--- $B9TNs7W;;%i%$%V%i%j$N=i4|2=(B. 
  !--- $B0z?t$O%9%+%i!<$N(B z $BJ}8~$N7W;;NN0h$NBg$-$5(B
  call linlib_init(cfgfile, RegZMax - RegZMin)
  

!!!
!!!$BJQ?t$N=i4|2=(B
!!!  
  !--- $BG[Ns$N3d$jEv$F(B, $B=i4|2=(B
  call ArareAlloc

  !--- $B>qMp>l$N=i4|2=(B
  ss_Exner_nl   = ss_ExnerDisturb
  ss_PotTemp_nl = ss_PotTempDisturb
  fs_VelX_nl    = fs_VelXDisturb
  sf_VelZ_nl    = sf_VelZDisturb
  
  !=== $B=i4|$N2aK0OBEY$r7W;;(B
  call SuperSaturationRatio(ss_Exner_nl,  &      !(in)
    &                       ss_PotTemp_nl, &     !(in)
    &                       ss_SupSatRatio_nl)   !(out)

!!!
!!!$B?tCM@QJ,(B. 1 $B%9%F%C%WJ,$@$12s$9(B. 
!!!
  !=== $BC10L6E7kNL$"$?$j$N@xG.$N7W;;(B
  call LatentHeatPerMass(ss_LatHeatPerMass_nl)   !(out)

  !--- $B12G4@-78?t(B, $B123H;678?t$N7W;;(B
  call EddyViscosity(                                     &
    & ss_Km_nl,                                        & !$B=i4|CM$O%<%m(B
    & DelTimeLong,                                     & 
    & fs_VelX_nl, sf_VelZ_nl, ss_Km_nl,                & !$B0\N.(B
    & fs_VelX_nl, sf_VelZ_nl, ss_PotTemp_nl, ss_Km_nl, & !$B@8@.(B, $B3H;6(B, $B>C;6(B
    & ss_Km_al, ss_Kh_al )

  !--- $BD9$$;~4V%9%F%C%W$G$N290L$N7W;;(B
  call Scalar(                                  & 
    & ss_PotTemp_nl,                         &
    & DelTimeLong,                           &
    & fs_VelX_nl, sf_VelZ_nl, ss_PotTemp_nl, & !$B0\N.(B
    & ss_PotTemp_nl, ss_Kh_nl,               & !$B@8@.(B, $BMpN.3H;6(B, $B?tCMG4@-(B
    & ss_PotTemp_al,                         & !(out) $B<!$N;~4V$N290L(B
    & ss_PotTempTendAdv_nl,                  & !(out) $B0\N.9`(B
    & ss_PotTempTendSrc_nl,                  & !(out) $B@8@.9`(B
    & ss_PotTempTendTurb_nl,                 & !(out) $BMpN.3H;69`(B
    & ss_PotTempTendNumDiff_nl)                !(out) $B?tCM3H;69`(B
  
  !--- $BD9$$;~4V%9%F%C%W$G$NB.EY(B u $B$N0\N.(B, $B3H;6(B, $BG4@-(B
!  call AdvectDiffuseX(                     &
!    & fs_VelX_nl, sf_VelZ_nl,           & !$B0\N.(B
!    & fs_VelX_nl, sf_VelZ_nl, ss_Km_nl, & !$BMpN.3H;6(B, $B?tCMG4@-(B
!    & fs_AdvDiffX_nl,  &                  !$B0\N.!&3H;6$K$h$kB.EYJQ2=N((B
!    & fs_VelXTendAdv_nl,  &                  !$B0\N.$K$h$kB.EYJQ2=N((B
!    & fs_VelXTendTurb_nl, &                  !$BMpN.3H;6$K$h$kB.EYJQ2=N((B
!    & fs_VelXTendNumDiff_nl)                 !$B?tCM3H;6$K$h$kB.EYJQ2=N((B

  !--- $BD9$$;~4V%9%F%C%W$G$NB.EY(B w $B$N0\N.(B, $B3H;6(B, $BG4@-(B
!  call AdvectDiffuseZ(                          &
!    & fs_VelX_nl, sf_VelZ_nl, ss_PotTemp_nl, & !$B0\N.(B, $BIbNO(B
!    & fs_VelX_nl, sf_VelZ_nl, ss_Km_nl,      & !$BMpN.3H;6(B, $B?tCMG4@-(B
!    & sf_AdvDiffZ_nl,  &
!    & sf_VelZTendAdv_nl, sf_VelZTendBuoy_nl, & 
!    & sf_VelZTendTurb_nl, sf_VelZTendNumDiff_nl)
  
  !--- $BC;$$%?%$%`%9%F%C%W$N=i4|CM$r:n$k(B

  fs_VelX_ns  = fs_VelX_nl
  sf_VelZ_ns  = sf_VelZ_nl
  ss_Exner_ns = ss_Exner_nl 
  ss_PotTemp_ns = ss_PotTemp_nl 
  ss_DensCloud_ns = ss_DensCloud_nl 
  ss_SupSatRatio_ns = ss_SupSatRatio_nl 
  
  !--- $BC;$$%?%$%`%9%F%C%W$G$N@QJ,(B
  do tau = 1, NstepShort / 2
    !=== $B6E7kNL$N7W;;(B
!    call MassCondense(  &
!      & ss_LatHeatPerMass_nl, &  !(in)
!      & ss_SupSatRatio_ns, &     !(in)
!      & ss_PotTemp_ns,     &     !(in)
!      & ss_Exner_ns, &           !(in)
!      & ss_DensCloud_ns, &       !(in)
!      & ss_MassCond_ns)          !(out)
 
    !=== $B6E7kJ*<AG;EY$N7W;;(B
!    call DensityCloud(  &
!      & ss_DensCloud_ns,                          &
!      & DelTimeShort,                                & !$B;~4V4V3V(B
!      & fs_VelX_nl, sf_VelZ_nl, ss_DensCloud_nl,  & !$B0\N.7W;;$GMQ$$$kCM(B
!      & ss_MassCond_ns,                         & !$B6E7kNL(B
!      & ss_DensCloud_as)                          

    !--- $B290L$N7W;;(B
    !==== $B6E7k@xG.$N7W;;(B
!    call LatentHeat(ss_MassCond_ns, ss_LatHeatPerMass_nl, ss_QCond)
    
    !==== 1 $B%9%F%C%W@h$N290L$r7W;;(B
    ss_PotTemp_as = ss_PotTemp_ns  &
      & + ( &
      &     - ss_PotTempTendAdv_nl & 
      &   ) * DelTimeShort
    call boundary(ss_BC, ss_PotTemp_as)

    !--- $BB.EY(B u $B$N7W;;(B
!    call VelocityX(                                             &
!      & fs_VelX_ns,                                          &
!      & DelTimeShort,                                        &
!      & fs_AdvDiffX_nl, fs_VelX_ns, sf_VelZ_ns, ss_Exner_ns, &
!      & fs_VelX_as )
     fs_VelX_as = fs_VelX_ns
    
    !--- $B%(%/%9%J!<4X?t$N7W;;(B
    !==== $B6E7k$K$h$k05NOJQ2=$r7W;;(B
!    call DExnerDTimeCondense(   &
!      & ss_MassCond_ns,         &  !(in)  $B6E7kNL(B
!      & ss_LatHeatPerMass_nl,   &  !(in)  $BC10L<ANL$"$?$j$N@xG.(B
!      & ss_ExnerTendCond_ns)       !(out) $B6E7k$KH<$&05NOJQ2=(B

    !==== $B05NO$N7W;;(B
!    call Exner( &
!      & sf_AdvDiffZ_nl, ss_ExnerTendCond_ns,  &              !(in)
!      & fs_VelX_ns, fs_VelX_as, sf_VelZ_ns, ss_Exner_ns, &   !(in)
!      & ss_Exner_as )                                        !(out)
     ss_Exner_as = ss_Exner_ns
    
    !--- $BB.EY(B w $B$N7W;;(B
!    call VelocityZ(  &
!      & sf_VelZ_nl, &
!      & DelTimeShort,    &
!      & sf_AdvDiffZ_nl, fs_VelX_ns, sf_VelZ_ns, ss_Exner_ns, ss_Exner_as, &
!      & sf_VelZ_as,  &
!      & sf_VelZTendPressGrad_ns)
     sf_VelZ_as = sf_VelZ_ns

    !=== $B2aK0OBEY$r7W;;(B
    call SuperSaturationRatio(ss_Exner_as,  &      !(in)
      &                       ss_PotTemp_as, &     !(in)
      &                       ss_SupSatRatio_as)   !(out)

    ss_Exner_as  = ss_Exner_ns 
    fs_VelX_as   = fs_VelX_ns
    sf_VelZ_as   = sf_VelZ_ns

    !--- $B%k!<%W$r2s$9$?$a$N=hCV(B
    ss_Exner_ns      = ss_Exner_as 
    ss_PotTemp_ns    = ss_PotTemp_as 
    ss_DensCloud_ns  = ss_DensCloud_as 
    fs_VelX_ns       = fs_VelX_as
    sf_VelZ_ns       = sf_VelZ_as
    ss_SupSatRatio_ns  = ss_SupSatRatio_as 

  end do
     
  !--- $B:G=*E*$JB.EY$rD9$$;~4V%9%F%C%W$G$N$b$N$H$_$J$9(B
  ss_Exner_al      = ss_Exner_as 
  ss_PotTemp_al    = ss_PotTemp_as 
  ss_DensCloud_al  = ss_DensCloud_as 
  fs_VelX_al       = fs_VelX_as
  sf_VelZ_al       = sf_VelZ_as     
  ss_SupSatRatio_al  = ss_SupSatRatio_as 

  !--- $B%U%!%$%k=PNO(B
  ReStartTime(1) = 0.0d0
  ReStartTime(2) = DelTimeLong
  call MakeReStartFile( ReStartFile, ReStartTime,  &
    & fs_VelX_nl,    sf_VelZ_nl,  ss_Exner_nl,  &
    & ss_PotTemp_nl, ss_DensCloud_nl,   &
    & ss_Km_nl,      ss_Kh_nl,    ss_SupSatRatio_nl, &
    & fs_VelX_al,    sf_VelZ_al,  ss_Exner_al,  &
    & ss_PotTemp_al, ss_DensCloud_al,  & 
    & ss_Km_al,      ss_Kh_al,    ss_SupSatRatio_al) 

contains


!!!
!!!$BG[Ns$N3d$jEv$F(B, $B=i4|2=(B
!!!
  subroutine ArareAlloc

    call BeginSub("ArareAlloc", &
      &            fmt="%c",     &
      &            c1="Allocate and initialize allocatable variables.")

    allocate(ss_Exner_as(DimXMin:DimXMax, DimZMin:DimZMax),  &
      & ss_PotTemp_as(DimXMin:DimXMax, DimZMin:DimZMax),  &
      & ss_DensCloud_as(DimXMin:DimXMax, DimZMin:DimZMax),  &
      & fs_VelX_as(DimXMin:DimXMax, DimZMin:DimZMax)  , &
      & sf_VelZ_as(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Exner_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_PotTemp_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_DensCloud_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelX_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZ_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelX_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZ_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Exner_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_PotTemp_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelX_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZ_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Exner_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_PotTemp_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_AdvDiffX_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_AdvDiffZ_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_ExnerTendCond_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelXTendAdv_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelXTendTurb_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & fs_VelXTendNumDiff_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZTendPressGrad_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZTendAdv_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZTendTurb_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZTendBuoy_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & sf_VelZTendNumDiff_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Km_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Kh_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Km_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_Kh_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_LatHeatPerMass_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_MassCond_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_DensCloud_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_DensCloud_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_SupSatRatio_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_SupSatRatio_al(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_SupSatRatio_ns(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_SupSatRatio_as(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_QCond(DimXMin:DimXMax, DimZMin:DimZMax))

    allocate( ss_PotTempTendAdv_nl(DimXMin:DimXMax, DimZMin:DimZMax),  &
      & ss_PotTempTendSrc_nl(DimXMin:DimXMax, DimZMin:DimZMax),  &
      & ss_PotTempTendTurb_nl(DimXMin:DimXMax, DimZMin:DimZMax), &
      & ss_PotTempTendNumDiff_nl(DimXMin:DimXMax, DimZMin:DimZMax))

    fs_VelX_al = 0.0d0
    fs_VelX_nl = 0.0d0

    sf_VelZ_al = 0.0d0
    sf_VelZ_nl = 0.0d0

    ss_Exner_al = 0.0d0
    ss_Exner_nl = 0.0d0

    ss_PotTemp_al = 0.0d0
    ss_PotTemp_nl = 0.0d0

    fs_VelX_as = 0.0d0
    fs_VelX_ns = 0.0d0

    sf_VelZ_as = 0.0d0
    sf_VelZ_ns = 0.0d0

    ss_Exner_as = 0.0d0
    ss_Exner_ns = 0.0d0

    ss_PotTemp_as = 0.0d0
    ss_PotTemp_ns = 0.0d0

    ss_DensCloud_as = 0.0d0
    ss_DensCloud_ns = 0.0d0

    fs_AdvDiffX_nl = 0.0d0
    sf_AdvDiffZ_nl = 0.0d0
    ss_ExnerTendCond_ns = 0.0d0

    ss_Km_nl = 0.0d0
    ss_Kh_nl = 0.0d0

    ss_Km_al = 0.0d0
    ss_Kh_al = 0.0d0

    ss_LatHeatPerMass_nl = 0.0d0
    ss_MassCond_ns = 0.0d0

    ss_DensCloud_nl = 0.0d0
    ss_DensCloud_al = 0.0d0

    ss_SupSatRatio_nl = 0.0d0
    ss_SupSatRatio_al = 0.0d0
    ss_SupSatRatio_ns = 0.0d0
    ss_SupSatRatio_as = 0.0d0

    ss_PotTempTendAdv_nl  = 0.0d0
    ss_PotTempTendSrc_nl  = 0.0d0
    ss_PotTempTendTurb_nl = 0.0d0
    ss_PotTempTendNumDiff_nl = 0.0d0

    ss_QCond = 0.0d0

    call EndSub("ArareAlloc")

  end subroutine ArareAlloc


end program ArareInit
