xref: /phasta/phSolver/common/bctint.f (revision 10d9781b35c469ead1916fa38f8accccbb97d1c2)
159599516SKenneth E. Jansenc-----------------------------------------------------------------------
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc  This module conveys temporal BC data.  Below functions read in the data
459599516SKenneth E. Jansenc  and interpolate it to the current time level.
559599516SKenneth E. Jansenc
659599516SKenneth E. Jansenc-----------------------------------------------------------------------
759599516SKenneth E. Jansen      module specialBC
859599516SKenneth E. Jansen      use pointer_data
959599516SKenneth E. Jansen
1059599516SKenneth E. Jansen      real*8, allocatable ::  BCt(:,:), acs(:,:), spamp(:)
1159599516SKenneth E. Jansen      real*8, allocatable ::  ytarget(:,:)
1259599516SKenneth E. Jansen      real*8, allocatable ::  BCperiod(:)
1359599516SKenneth E. Jansen      integer, allocatable :: nBCt(:), numBCt(:)
1459599516SKenneth E. Jansen
1559599516SKenneth E. Jansen      type(r2d),dimension(:),allocatable :: BCtptr
1659599516SKenneth E. Jansen        ! BCtptr%p(:,:) is to replace BCt(ic,:,:), in which
1759599516SKenneth E. Jansen        ! the second index is dynamically decided in
1859599516SKenneth E. Jansen        ! the setup stage.
1959599516SKenneth E. Jansen
2059599516SKenneth E. Jansen      integer ntv,nptsmax
2159599516SKenneth E. Jansenc$$$      integer itvn
2259599516SKenneth E. Jansen      end module
2359599516SKenneth E. Jansen
2459599516SKenneth E. Jansenc-----------------------------------------------------------------------
2559599516SKenneth E. Jansenc
2659599516SKenneth E. Jansenc  This module conveys flow rate history for the different impedance outlets
2759599516SKenneth E. Jansenc  over one period. Below functions read in the data and store it for the
2859599516SKenneth E. Jansenc  current time level.
2959599516SKenneth E. Jansenc
3059599516SKenneth E. Jansenc-----------------------------------------------------------------------
3159599516SKenneth E. Jansen      module convolImpFlow
3259599516SKenneth E. Jansen
3359599516SKenneth E. Jansen      real*8, allocatable ::  QHistImp(:,:), ValueImpt(:,:,:)
3459599516SKenneth E. Jansen      real*8, allocatable ::  ValueListImp(:,:), ConvCoef(:,:)
3559599516SKenneth E. Jansen      real*8, allocatable ::  ImpConvCoef(:,:), poldImp(:)
3659599516SKenneth E. Jansen      integer ntimeptpT,numDataImp
3759599516SKenneth E. Jansen      integer, allocatable :: nImpt(:), numImpt(:)
3859599516SKenneth E. Jansen      integer nptsImpmax
3959599516SKenneth E. Jansen      real*8, allocatable ::  QHistTry(:,:), QHistTryF(:,:) !filter
4059599516SKenneth E. Jansen      integer cutfreq !filter
4159599516SKenneth E. Jansen      end module
4259599516SKenneth E. Jansenc-----------------------------------------------------------------------
4359599516SKenneth E. Jansenc
4459599516SKenneth E. Jansenc  This module conveys the parameters for the different RCR outlets.
4559599516SKenneth E. Jansenc  Below functions read in the inputs (proximal resistance, capacitance,
4659599516SKenneth E. Jansenc  distal resistance and distal pressure) and store it for the
4759599516SKenneth E. Jansenc  current time level.
4859599516SKenneth E. Jansenc
4959599516SKenneth E. Jansenc-----------------------------------------------------------------------
5059599516SKenneth E. Jansen      module convolRCRFlow
5159599516SKenneth E. Jansen
5259599516SKenneth E. Jansen      real*8, allocatable ::  ValueListRCR(:,:), ValuePdist(:,:,:) !inputs
5359599516SKenneth E. Jansen      real*8, allocatable ::  QHistRCR(:,:), HopRCR(:) !calc
5459599516SKenneth E. Jansen      real*8, allocatable ::  RCRConvCoef(:,:), poldRCR(:) !calc
5559599516SKenneth E. Jansen      real*8, allocatable ::  dtRCR(:) !scaled timestep: deltat/RdC
5659599516SKenneth E. Jansen      real*8, allocatable ::  RCRic(:) !(P(0)-RQ(0)-Pd(0))
5759599516SKenneth E. Jansen      integer nptsRCRmax,numDataRCR !to read inputs
5859599516SKenneth E. Jansen      integer, allocatable :: numRCRt(:) !to read inputs
5959599516SKenneth E. Jansen      end module
6059599516SKenneth E. Jansen
6159599516SKenneth E. Jansenc-----------------------------------------------------------------------
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansenc     Initialize:
6459599516SKenneth E. Jansenc
6559599516SKenneth E. Jansenc-----------------------------------------------------------------------
6659599516SKenneth E. Jansen      subroutine initSponge( y,x)
6759599516SKenneth E. Jansen
6859599516SKenneth E. Jansen      use     specialBC
6959599516SKenneth E. Jansen      include "common.h"
7059599516SKenneth E. Jansen
7159599516SKenneth E. Jansen      real*8   y(nshg,nflow), x(numnp,3)
7259599516SKenneth E. Jansen      allocate (ytarget(nshg,nflow))
7359599516SKenneth E. Jansen
7459599516SKenneth E. Jansen      if(matflg(5,1).eq.5) then
7559599516SKenneth E. Jansen         write(*,*) 'calculating IC sponge'
7659599516SKenneth E. Jansen         ytarget = y
7759599516SKenneth E. Jansen      else
78*10d9781bSKenneth E. Jansen         ytarget(:,1)= xvel_ini
79*10d9781bSKenneth E. Jansen         ytarget(:,2)= yvel_ini
80*10d9781bSKenneth E. Jansen         ytarget(:,3)= zvel_ini
81*10d9781bSKenneth E. Jansen         ytarget(:,4)= pres_ini
82*10d9781bSKenneth E. Jansen         ytarget(:,5)= temp_ini
83*10d9781bSKenneth E. Jansen         if(myrank.eq.0)write(*,*) 'calculating Analytic sponge'
84*10d9781bSKenneth E. Jansen       endif
85*10d9781bSKenneth E. Jansen!      else
86*10d9781bSKenneth E. Jansen       if (0.eq.1) then
8759599516SKenneth E. Jansen
8859599516SKenneth E. Jansenc
8959599516SKenneth E. Jansenc OLD style sponge pushed onto target.  You need to be sure that your
9059599516SKenneth E. Jansenc solver.inp entries for start and stop of sponge match as well as the
9159599516SKenneth E. Jansenc growth rates
9259599516SKenneth E. Jansenc
9359599516SKenneth E. Jansen      vcl=datmat(1,5,1)         ! velocity on centerline
9459599516SKenneth E. Jansen      rslc=datmat(2,5,1)        ! shear layer center radius
9559599516SKenneth E. Jansen      bfz=datmat(3,5,1)
9659599516SKenneth E. Jansen      we=3.0*29./682.
9759599516SKenneth E. Jansen      rsteep=3.0
9859599516SKenneth E. Jansen      zstart=30.0
9959599516SKenneth E. Jansen      radst=10.0
10059599516SKenneth E. Jansen      radsts=radst*radst
10159599516SKenneth E. Jansen      do id=1,numnp
10259599516SKenneth E. Jansen         radsqr=x(id,2)**2+x(id,1)**2
10359599516SKenneth E. Jansenc         if((x(id,3).gt. zstart) .or. (radsqr.gt.radsts))  then
10459599516SKenneth E. Jansen            rad=sqrt(radsqr)
10559599516SKenneth E. Jansen            radc=max(rad,radst)
10659599516SKenneth E. Jansen            zval=max(x(id,3),zstart)
10759599516SKenneth E. Jansen            utarget=(tanh(rsteep*(rslc-rad))+one)/two*
10859599516SKenneth E. Jansen     &                    (vcl-we) + we
10959599516SKenneth E. Jansen            Ttarget  = press/(ro*Rgas)
11059599516SKenneth E. Jansen            ptarget= press
11159599516SKenneth E. Jansen            ytarget(id,1) = zero
11259599516SKenneth E. Jansen            ytarget(id,2) = zero
11359599516SKenneth E. Jansen            ytarget(id,3) = utarget
11459599516SKenneth E. Jansen            ytarget(id,4) = ptarget
11559599516SKenneth E. Jansen            ytarget(id,5) = Ttarget
11659599516SKenneth E. Jansenc         endif
11759599516SKenneth E. Jansen      enddo
11859599516SKenneth E. Jansen      endif
11959599516SKenneth E. Jansen      return
12059599516SKenneth E. Jansen      end
12159599516SKenneth E. Jansen
12259599516SKenneth E. Jansen
12359599516SKenneth E. Jansenc-----------------------------------------------------------------------
12459599516SKenneth E. Jansenc
12559599516SKenneth E. Jansenc     Initialize:time varying boundary condition
12659599516SKenneth E. Jansenc
12759599516SKenneth E. Jansenc-----------------------------------------------------------------------
12859599516SKenneth E. Jansen      subroutine initBCt( x, iBC, BC )
12959599516SKenneth E. Jansen
13059599516SKenneth E. Jansen      use     specialBC
13159599516SKenneth E. Jansen      include "common.h"
13259599516SKenneth E. Jansen
13359599516SKenneth E. Jansen      real*8   x(numnp,nsd), BC(nshg,ndofBC), rj1,rj2,rj3,rj4,distd,epsd
13459599516SKenneth E. Jansen      integer  iBC(nshg)
13559599516SKenneth E. Jansen      character*80 card
13659599516SKenneth E. Jansen      real*8 distds
13759599516SKenneth E. Jansen      real*8 dd
13859599516SKenneth E. Jansen
13959599516SKenneth E. Jansen      integer irstart
14059599516SKenneth E. Jansen      real(kind=8),allocatable,dimension(:) ::  bcttest
14159599516SKenneth E. Jansen
14259599516SKenneth E. Jansen      real*8 t0,tlen,t1,tstart,tend
14359599516SKenneth E. Jansen      integer i0,ilen,i1,nper,istart,iend
14459599516SKenneth E. Jansenc
14559599516SKenneth E. Jansenc  This one should be used for boundary layer meshes where bct.dat must
14659599516SKenneth E. Jansenc  be given to greater precision than is currently being generated.
14759599516SKenneth E. Jansenc
14859599516SKenneth E. Jansenc      epsd=1.0d-12    ! this is distance SQUARED to save square root
14959599516SKenneth E. Jansen
15059599516SKenneth E. Jansen      epsd=1.0d-8              ! this is distance SQUARED to save square root
15159599516SKenneth E. Jansen
15259599516SKenneth E. Jansen      ic=0                      !count the number on this processor
15359599516SKenneth E. Jansen
15459599516SKenneth E. Jansen! ************ Chun Sun: limit the memory use if required time steps
15559599516SKenneth E. Jansen! ************ is smaller than total bct.dat length.
15659599516SKenneth E. Jansen
15759599516SKenneth E. Jansen      ! readin the start timestep
15859599516SKenneth E. Jansen
15959599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
16059599516SKenneth E. Jansen      read(72,*) irstart
16159599516SKenneth E. Jansen      close(72)
16259599516SKenneth E. Jansen
16359599516SKenneth E. Jansen      ! use irstart-1 and nstep+1 to avoid machine tolerance issues
16459599516SKenneth E. Jansen      t0 = max(zero,(irstart-1)*Delt(1))
16559599516SKenneth E. Jansen      tlen = (nstep(1)+1)*Delt(1)
16659599516SKenneth E. Jansen      ! Assumption: constant time step in one run. If you use variable
16759599516SKenneth E. Jansen      ! time step in one run, this should be modified.
16859599516SKenneth E. Jansen      t1 = t0 + tlen
16959599516SKenneth E. Jansen
17059599516SKenneth E. Jansen      if(myrank.eq.master)
17159599516SKenneth E. Jansen     &   write(*,*) 'necessary bct timing: from ',t0,' to ',t1
17259599516SKenneth E. Jansen
17359599516SKenneth E. Jansen! **************************************************************
17459599516SKenneth E. Jansen
17559599516SKenneth E. Jansen
17659599516SKenneth E. Jansen      if(any(ibits(iBC,3,3).eq.7)) then
17759599516SKenneth E. Jansen         if(myrank.eq.master) write(*,*) 'opening bct.dat'
17859599516SKenneth E. Jansenc         open(unit=567, file='bct.dat',status='old')
17959599516SKenneth E. Jansen         open(unit=567, file='bct.dat',ACTION='READ',STATUS='old')
18059599516SKenneth E. Jansen         read (567,'(a80)') card
18159599516SKenneth E. Jansen           read (card,*) ntv, nptsmax
18259599516SKenneth E. Jansenc        read(567,*) ntv,nptsmax
18359599516SKenneth E. Jansen         allocate (nBCt(numnp))
18459599516SKenneth E. Jansen         allocate (numBCt(ntv))
18559599516SKenneth E. Jansen         allocate (BCt(nptsmax,4))
18659599516SKenneth E. Jansen         allocate (BCperiod(ntv))
18759599516SKenneth E. Jansen         allocate (BCtptr(ntv))  ! dynamic bct allocation
18859599516SKenneth E. Jansen         do k=1,ntv
18959599516SKenneth E. Jansen            read(567,*) x1,x2,x3,ntpts
19059599516SKenneth E. Jansenc
19159599516SKenneth E. Jansenc Find the point on the boundary (if it is on this processor)
19259599516SKenneth E. Jansenc that matches this point
19359599516SKenneth E. Jansenc
19459599516SKenneth E. Jansen            do i=1,numnp
19559599516SKenneth E. Jansen               if(ibits(ibc(i),3,3) .eq.7) then
19659599516SKenneth E. Jansen                  dd= distds(x1,x2,x3,x(i,1),x(i,2),x(i,3))
19759599516SKenneth E. Jansen                  if(dd.lt.epsd) then
19859599516SKenneth E. Jansen                     ic=ic+1
19959599516SKenneth E. Jansen                     nBCt(ic)=i ! the pointer to this point
20059599516SKenneth E. Jansen!                     numBCt(ic)=ntpts ! the number of time series
20159599516SKenneth E. Jansen                     do j=1,ntpts
20259599516SKenneth E. Jansenc                        read(567,*) BCt(ic,j,4),(BCt(ic,j,n),n=1,3)
20359599516SKenneth E. Jansen                        read(567,*) (BCt(j,n),n=1,4)
20459599516SKenneth E. Jansen                     enddo
20559599516SKenneth E. Jansen              ! at this point we have all bc data of spacial point
20659599516SKenneth E. Jansen              ! ic/i in all time domain. now we figure out which subset
20759599516SKenneth E. Jansen              ! is necessary.
20859599516SKenneth E. Jansen                     if (tlen.ge.BCt(ntpts,4)) then
20959599516SKenneth E. Jansen                         ! whole run is larger than whole period
21059599516SKenneth E. Jansen                         ! should take all data
21159599516SKenneth E. Jansen                         ilen = ntpts
21259599516SKenneth E. Jansen                         allocate(BCtptr(ic)%p(ilen,4))
21359599516SKenneth E. Jansen                         BCtptr(ic)%p = BCt(1:ilen,:)
21459599516SKenneth E. Jansen                     else
21559599516SKenneth E. Jansen                         nper = t0/BCt(ntpts,4)
21659599516SKenneth E. Jansen                         tstart = t0-nper*BCt(ntpts,4)
21759599516SKenneth E. Jansen                         nper = t1/BCt(ntpts,4)
21859599516SKenneth E. Jansen                         tend = t1-nper*BCt(ntpts,4)
21959599516SKenneth E. Jansen                         istart = iBfind(BCt(:,4),ntpts,tstart)
22059599516SKenneth E. Jansen                         iend = iBfind(BCt(:,4),ntpts,tend)
22159599516SKenneth E. Jansen                         if (istart>1) istart=istart-1
22259599516SKenneth E. Jansen                         if (iend<ntpts) iend=iend+1
22359599516SKenneth E. Jansen                         !write(*,*)'bcstart:',BCt(istart,4),tstart
22459599516SKenneth E. Jansen                         !write(*,*)'bcend:',BCt(iend,4),tend
22559599516SKenneth E. Jansen                         if (tstart.lt.tend) then ! does not loop
22659599516SKenneth E. Jansen                             ilen = iend-istart+1
22759599516SKenneth E. Jansen                             allocate(BCtptr(ic)%p(ilen,4))
22859599516SKenneth E. Jansen                             BCtptr(ic)%p = BCt(istart:iend,:)
22959599516SKenneth E. Jansen                         else ! loop within two consequetive time period
23059599516SKenneth E. Jansen                             i0 = (ntpts-istart+1) ! first segment
23159599516SKenneth E. Jansen                             ilen = iend + i0
23259599516SKenneth E. Jansen                             allocate(BCtptr(ic)%p(ilen,4))
23359599516SKenneth E. Jansen                      BCtptr(ic)%p(1:i0,:) = BCt(istart:ntpts,:)
23459599516SKenneth E. Jansen                      BCtptr(ic)%p(i0+1:ilen,:) = BCt(1:iend,:)
23559599516SKenneth E. Jansen                      BCtptr(ic)%p(i0+1:ilen,4) =
23659599516SKenneth E. Jansen     &                 BCtptr(ic)%p(i0+1:ilen,4) + BCt(ntpts,4)
23759599516SKenneth E. Jansen                         endif
23859599516SKenneth E. Jansen                     endif
23959599516SKenneth E. Jansen                     numBCt(ic) = ilen
24059599516SKenneth E. Jansen              BCtptr(ic)%p(:,4) = BCtptr(ic)%p(:,4)*bcttimescale
24159599516SKenneth E. Jansen              BCperiod(ic) = BCt(ntpts,4)*bcttimescale
24259599516SKenneth E. Jansen                     exit
24359599516SKenneth E. Jansen                  endif
24459599516SKenneth E. Jansen               endif
24559599516SKenneth E. Jansen            enddo
24659599516SKenneth E. Jansen            if(i.eq.numnp+1) then
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansenc  if we get here the point was not found.  It must be on another
24959599516SKenneth E. Jansenc  processor so we read past this record and move on
25059599516SKenneth E. Jansenc
25159599516SKenneth E. Jansen               do j=1,ntpts
25259599516SKenneth E. Jansen                  read(567,*) rj1,rj2,rj3,rj4
25359599516SKenneth E. Jansen               enddo
25459599516SKenneth E. Jansen            endif
25559599516SKenneth E. Jansen         enddo                  ! end of the loop over ntv
25659599516SKenneth E. Jansen         !BCt(:,:,4)=BCt(:,:,4)*bcttimescale
25759599516SKenneth E. Jansen         deallocate(BCt)
25859599516SKenneth E. Jansen      endif                     ! any 3 component nodes
25959599516SKenneth E. Jansen      itvn=ic
26059599516SKenneth E. Jansen      close(567)
26159599516SKenneth E. Jansen      if (ic.gt.0) then
26259599516SKenneth E. Jansen         write(*,*)'myrank=',myrank,' and I found ',ic,' nodes.'
26359599516SKenneth E. Jansenc      else
26459599516SKenneth E. Jansenc         deallocate(nBCt)
26559599516SKenneth E. Jansenc         deallocate(numBCt)
26659599516SKenneth E. Jansenc         deallocate(BCt)
26759599516SKenneth E. Jansen      endif
26859599516SKenneth E. Jansen
26959599516SKenneth E. Jansen      return
27059599516SKenneth E. Jansen      end
27159599516SKenneth E. Jansen
27259599516SKenneth E. Jansen!***************************************************************
27359599516SKenneth E. Jansen!      A Binary search return an index of a real array.
27459599516SKenneth E. Jansen!      This array should be an ascending sorted array.
27559599516SKenneth E. Jansen!***************************************************************
27659599516SKenneth E. Jansen      function iBfind(bArray,bLen,bVal)
27759599516SKenneth E. Jansen      integer,intent(in) :: bLen
27859599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(bLen) :: bArray
27959599516SKenneth E. Jansen      real(kind=8),intent(in) :: bVal
28059599516SKenneth E. Jansen      integer mlen,indx,bstart,bend
28159599516SKenneth E. Jansen      bstart = 1
28259599516SKenneth E. Jansen      bend = bLen
28359599516SKenneth E. Jansen      indx = (bLen+1)/2
28459599516SKenneth E. Jansen      do while ((bstart+1).lt.bend)
28559599516SKenneth E. Jansen         if (bVal.gt.bArray(indx)) then
28659599516SKenneth E. Jansen             bstart = indx
28759599516SKenneth E. Jansen         else
28859599516SKenneth E. Jansen             bend = indx
28959599516SKenneth E. Jansen         endif
29059599516SKenneth E. Jansen         indx = (bstart+bend)/2
29159599516SKenneth E. Jansen      enddo
29259599516SKenneth E. Jansen      iBfind = indx
29359599516SKenneth E. Jansen      return
29459599516SKenneth E. Jansen      end function
29559599516SKenneth E. Jansen
29659599516SKenneth E. Jansen
29759599516SKenneth E. Jansen      subroutine BCint(timel,shp,shgl,shpb,shglb,x,BC,iBC)
29859599516SKenneth E. Jansen
29959599516SKenneth E. Jansen      use     specialBC ! brings in itvn,nbct, bct, numbct, nptsmax
30059599516SKenneth E. Jansen
30159599516SKenneth E. Jansen      include "common.h"
30259599516SKenneth E. Jansen
30359599516SKenneth E. Jansen      real*8   BC(nshg,ndofBC), timel,t
30459599516SKenneth E. Jansen      real*8   x(numnp,nsd),
30559599516SKenneth E. Jansen     &         shp(MAXTOP,maxsh,MAXQPT),
30659599516SKenneth E. Jansen     &         shgl(MAXTOP,nsd,maxsh,MAXQPT),
30759599516SKenneth E. Jansen     &         shpb(MAXTOP,maxsh,MAXQPT),
30859599516SKenneth E. Jansen     &         shglb(MAXTOP,nsd,maxsh,MAXQPT)
30959599516SKenneth E. Jansen
31059599516SKenneth E. Jansen      integer  iBC(numnp),nlast,i,j,nper
31159599516SKenneth E. Jansen
31259599516SKenneth E. Jansen      do i =1,itvn ! itvn is the number of varying nodes on this proc
31359599516SKenneth E. Jansen
31459599516SKenneth E. Jansen         nlast=numBCt(i)     ! number of time series to interpolate from
31559599516SKenneth E. Jansen         nper=timel/BCperiod(i)
31659599516SKenneth E. Jansen        ! number of periods completed to shift off
31759599516SKenneth E. Jansen
31859599516SKenneth E. Jansen
31959599516SKenneth E. Jansen         t=timel-nper*BCperiod(i)  ! now time in periodic domain
32059599516SKenneth E. Jansen
32159599516SKenneth E. Jansen         do j=2,nlast   !loop to find the interval that we are in
32259599516SKenneth E. Jansen
32359599516SKenneth E. Jansen            if(BCtptr(i)%p(j,4).gt.t) then
32459599516SKenneth E. Jansen            ! this is upper bound, j-1 is lower
32559599516SKenneth E. Jansen
32659599516SKenneth E. Jansen      wr=(t-BCtptr(i)%p(j-1,4))/(BCtptr(i)%p(j,4)-BCtptr(i)%p(j-1,4))
32759599516SKenneth E. Jansen               BC(nbct(i),3:5)= BCtptr(i)%p(j-1,1:3)*(one-wr)
32859599516SKenneth E. Jansen     &                        + BCtptr(i)%p(j,1:3)*wr
32959599516SKenneth E. Jansen               exit
33059599516SKenneth E. Jansen
33159599516SKenneth E. Jansen            endif
33259599516SKenneth E. Jansen         enddo
33359599516SKenneth E. Jansen      enddo
33459599516SKenneth E. Jansen      return
33559599516SKenneth E. Jansen      end
33659599516SKenneth E. Jansen
33759599516SKenneth E. Jansen      function distds(x1,y1,z1,x2,y2,z2)
33859599516SKenneth E. Jansen      real*8 distds
33959599516SKenneth E. Jansen      real*8 x1,y1,z1,x2,y2,z2,x,y,z
34059599516SKenneth E. Jansen      x=x1-x2
34159599516SKenneth E. Jansen      y=y1-y2
34259599516SKenneth E. Jansen      z=z1-z2
34359599516SKenneth E. Jansen      distds=x*x+y*y+z*z
34459599516SKenneth E. Jansen      return
34559599516SKenneth E. Jansen      end
34659599516SKenneth E. Jansenc-----------------------------------------------------------------------
34759599516SKenneth E. Jansenc   initialize the impedance boundary condition:
34859599516SKenneth E. Jansenc   read the data in initImpt
34959599516SKenneth E. Jansenc   interpolate the data to match the process time step in Impint
35059599516SKenneth E. Jansenc-----------------------------------------------------------------------
35159599516SKenneth E. Jansen      subroutine initImpt()
35259599516SKenneth E. Jansen
35359599516SKenneth E. Jansen      use convolImpFlow
35459599516SKenneth E. Jansen      include "common.h"
35559599516SKenneth E. Jansen
35659599516SKenneth E. Jansen      open(unit=817, file='impt.dat',status='old')
35759599516SKenneth E. Jansen         read (817,*) nptsImpmax
35859599516SKenneth E. Jansen         allocate (numImpt(numImpSrfs))
35959599516SKenneth E. Jansen         allocate (ValueImpt(nptsImpmax,2,numImpSrfs))
36059599516SKenneth E. Jansen         ValueImpt=0
36159599516SKenneth E. Jansen         do k=1,numImpSrfs
36259599516SKenneth E. Jansen            read (817,*) numDataImp
36359599516SKenneth E. Jansen            numImpt(k) = numDataImp
36459599516SKenneth E. Jansen            do j=1,numDataImp
36559599516SKenneth E. Jansen               read(817,*) (ValueImpt(j,n,k),n=1,2) ! n=1 time, 2 value
36659599516SKenneth E. Jansen            enddo
36759599516SKenneth E. Jansen         enddo
36859599516SKenneth E. Jansen      close(817)
36959599516SKenneth E. Jansen
37059599516SKenneth E. Jansen      allocate (ValueListImp(ntimeptpT+1,numImpSrfs))
37159599516SKenneth E. Jansen      ValueListImp = zero
37259599516SKenneth E. Jansen      ValueListImp(ntimeptpT+1,:) = ValueImpt(1,2,:) !Z(time=0), last entry
37359599516SKenneth E. Jansen      ValueListImp(1,:) = ValueImpt(1,2,:) !Z(time=0)=Z(time=T)
37459599516SKenneth E. Jansen      return
37559599516SKenneth E. Jansen      end
37659599516SKenneth E. Jansen
37759599516SKenneth E. Jansen
37859599516SKenneth E. Jansen
37959599516SKenneth E. Jansen      subroutine Impint(ctime,jstep)
38059599516SKenneth E. Jansen
38159599516SKenneth E. Jansen      use convolImpFlow
38259599516SKenneth E. Jansen      include "common.h"
38359599516SKenneth E. Jansen
38459599516SKenneth E. Jansen      real*8 ctime, ptime
38559599516SKenneth E. Jansen      integer nlast, nper, k, j , jstep
38659599516SKenneth E. Jansen
38759599516SKenneth E. Jansen
38859599516SKenneth E. Jansen      do k =1,numImpSrfs
38959599516SKenneth E. Jansen         nlast=numImpt(k)     ! number of time series to interpolate from
39059599516SKenneth E. Jansen         nper=ctime/ValueImpt(nlast,1,k)!number of periods completed to shift off
39159599516SKenneth E. Jansen         ptime = ctime-nper*ValueImpt(nlast,1,k)  ! now time in periodic domain
39259599516SKenneth E. Jansen
39359599516SKenneth E. Jansen         do j=2,nlast   !loop to find the interval that we are in
39459599516SKenneth E. Jansen
39559599516SKenneth E. Jansen            if(ValueImpt(j,1,k).gt.ptime) then  ! this is upper bound, j-1 is lower
39659599516SKenneth E. Jansen               wr=(ptime-ValueImpt(j-1,1,k))
39759599516SKenneth E. Jansen     &             / ( ValueImpt(j,1,k)-ValueImpt(j-1,1,k) )
39859599516SKenneth E. Jansen               ValueListImp(jstep,k)= ValueImpt(j-1,2,k)*(one-wr)
39959599516SKenneth E. Jansen     &                        + ValueImpt(j,2,k)*wr
40059599516SKenneth E. Jansen               exit
40159599516SKenneth E. Jansen            endif
40259599516SKenneth E. Jansen
40359599516SKenneth E. Jansen         enddo
40459599516SKenneth E. Jansen      enddo
40559599516SKenneth E. Jansen      return
40659599516SKenneth E. Jansen      end
40759599516SKenneth E. Jansen
40859599516SKenneth E. Jansenc-----------------------------------------------------------------------------
40959599516SKenneth E. Jansenc     time filter for a periodic function (sin cardinal + window function)
41059599516SKenneth E. Jansenc     is used for the impedance and the flow rate history
41159599516SKenneth E. Jansenc-----------------------------------------------------------------------------
41259599516SKenneth E. Jansen      subroutine Filter(Filtered,DataHist,nptf,timestep,cutfreq)
41359599516SKenneth E. Jansen
41459599516SKenneth E. Jansen      include "common.h"
41559599516SKenneth E. Jansen
41659599516SKenneth E. Jansen      integer nptf, cutfreq, j, k, m, s, Filtime(nptf)
41759599516SKenneth E. Jansen      real*8  DataHist(nptf,numImpSrfs), Window(nptf)
41859599516SKenneth E. Jansen      real*8  Sinc(nptf), FilterSW(nptf), Filtered(nptf,numImpSrfs)
41959599516SKenneth E. Jansen      real*8  windK, timestep
42059599516SKenneth E. Jansen
42159599516SKenneth E. Jansen      windK = cutfreq*2 + 1
42259599516SKenneth E. Jansen      do j=1,nptf
42359599516SKenneth E. Jansen         Filtime(j) = j-1
42459599516SKenneth E. Jansen         Window(j) = 0.42+0.5*cos(2*pi*Filtime(j)/nptf)
42559599516SKenneth E. Jansen     &              +0.08*cos(4*pi*Filtime(j)/nptf)
42659599516SKenneth E. Jansen         Sinc(j) = sin(pi*Filtime(j)*windK/nptf)/sin(pi*Filtime(j)/nptf)
42759599516SKenneth E. Jansen      enddo
42859599516SKenneth E. Jansen      Sinc(1) = windK
42959599516SKenneth E. Jansen
43059599516SKenneth E. Jansen      do j=1,nptf
43159599516SKenneth E. Jansen         FilterSW(j) = Window(nptf+1-j)*Sinc(nptf+1-j) !filter for convolution
43259599516SKenneth E. Jansen      enddo
43359599516SKenneth E. Jansen
43459599516SKenneth E. Jansen      Filtered = zero
43559599516SKenneth E. Jansen      do m=1,nptf
43659599516SKenneth E. Jansen         do j=1,nptf
43759599516SKenneth E. Jansen            s=modulo(m-nptf+j,nptf)
43859599516SKenneth E. Jansen            if(s.eq.zero) then
43959599516SKenneth E. Jansen               s=nptf
44059599516SKenneth E. Jansen            endif
44159599516SKenneth E. Jansen            Filtered(m,:) = Filtered(m,:)
44259599516SKenneth E. Jansen     &              +FilterSW(j)*DataHist(s,:)/nptf !filter convolution
44359599516SKenneth E. Jansen         enddo
44459599516SKenneth E. Jansen      enddo
44559599516SKenneth E. Jansen
44659599516SKenneth E. Jansen      return
44759599516SKenneth E. Jansen      end
44859599516SKenneth E. Jansenc-----------------------------------------------------------------------
44959599516SKenneth E. Jansenc   initialize the RCR boundary condition:
45059599516SKenneth E. Jansenc   read the data in initRCRt
45159599516SKenneth E. Jansenc   interpolate the data to match the process time step in RCRint
45259599516SKenneth E. Jansenc-----------------------------------------------------------------------
45359599516SKenneth E. Jansen      subroutine initRCRt()
45459599516SKenneth E. Jansen
45559599516SKenneth E. Jansen      use convolRCRFlow
45659599516SKenneth E. Jansen      include "common.h"
45759599516SKenneth E. Jansen
45859599516SKenneth E. Jansen      open(unit=818, file='rcrt.dat',status='old')
45959599516SKenneth E. Jansen         read (818,*) nptsRCRmax
46059599516SKenneth E. Jansen         allocate (numRCRt(numRCRSrfs))
46159599516SKenneth E. Jansen         allocate (ValuePdist(nptsRCRmax,2,numRCRSrfs))
46259599516SKenneth E. Jansen         allocate (ValueListRCR(3,numRCRSrfs))
46359599516SKenneth E. Jansen         ValuePdist=0
46459599516SKenneth E. Jansen         ValueListRCR=0
46559599516SKenneth E. Jansen         do k=1,numRCRSrfs
46659599516SKenneth E. Jansen            read (818,*) numDataRCR
46759599516SKenneth E. Jansen            numRCRt(k) = numDataRCR
46859599516SKenneth E. Jansen            do j=1,3
46959599516SKenneth E. Jansen               read(818,*) ValueListRCR(j,k) ! reads Rp,C,Rd
47059599516SKenneth E. Jansen            enddo
47159599516SKenneth E. Jansen            do j=1,numDataRCR
47259599516SKenneth E. Jansen               read(818,*) (ValuePdist(j,n,k),n=1,2) ! n=1 time, 2 value
47359599516SKenneth E. Jansen            enddo
47459599516SKenneth E. Jansen         enddo
47559599516SKenneth E. Jansen      close(818)
47659599516SKenneth E. Jansen
47759599516SKenneth E. Jansen      allocate (dtRCR(numRCRSrfs))
47859599516SKenneth E. Jansen
47959599516SKenneth E. Jansen      return
48059599516SKenneth E. Jansen      end
48159599516SKenneth E. Jansen
48259599516SKenneth E. Jansen
48359599516SKenneth E. Jansen      subroutine RCRint(ctime,Pdist)
48459599516SKenneth E. Jansen
48559599516SKenneth E. Jansen      use convolRCRFlow ! brings numRCRSrfs, ValuePdist
48659599516SKenneth E. Jansen      include "common.h"
48759599516SKenneth E. Jansen
48859599516SKenneth E. Jansen      real*8  ctime, ptime
48959599516SKenneth E. Jansen      integer nlast, nper, k, j
49059599516SKenneth E. Jansen      real*8  Pdist(0:MAXSURF)
49159599516SKenneth E. Jansen
49259599516SKenneth E. Jansen      do k =1,numRCRSrfs
49359599516SKenneth E. Jansen         nlast=numRCRt(k)     ! number of time series to interpolate from
49459599516SKenneth E. Jansen         nper=ctime/ValuePdist(nlast,1,k)!number of periods completed to shift off
49559599516SKenneth E. Jansen         ptime = ctime-nper*ValuePdist(nlast,1,k)  ! now time in periodic domain
49659599516SKenneth E. Jansen
49759599516SKenneth E. Jansen         do j=2,nlast   !loop to find the interval that we are in
49859599516SKenneth E. Jansen
49959599516SKenneth E. Jansen            if(ValuePdist(j,1,k).gt.ptime) then  ! this is upper bound, j-1 is lower
50059599516SKenneth E. Jansen               wr=(ptime-ValuePdist(j-1,1,k))
50159599516SKenneth E. Jansen     &             / ( ValuePdist(j,1,k)-ValuePdist(j-1,1,k) )
50259599516SKenneth E. Jansen               Pdist(k)= ValuePdist(j-1,2,k)*(one-wr)
50359599516SKenneth E. Jansen     &                        + ValuePdist(j,2,k)*wr
50459599516SKenneth E. Jansen               exit
50559599516SKenneth E. Jansen            endif
50659599516SKenneth E. Jansen
50759599516SKenneth E. Jansen         enddo
50859599516SKenneth E. Jansen      enddo
50959599516SKenneth E. Jansen      return
51059599516SKenneth E. Jansen      end
51159599516SKenneth E. Jansen
51259599516SKenneth E. Jansenc-----------------------------------------------------------------------
51359599516SKenneth E. Jansenc     returns in pold the history dependent part of the pressure in the
51459599516SKenneth E. Jansenc     impedance/flow rate convolution for the impedance and RCR BC
51559599516SKenneth E. Jansenc-----------------------------------------------------------------------
51659599516SKenneth E. Jansen      subroutine pHist(pressHist,QHist,betas,nTimePoint,nSrfs)
51759599516SKenneth E. Jansen
51859599516SKenneth E. Jansen      include "common.h"
51959599516SKenneth E. Jansen
52059599516SKenneth E. Jansen      integer  nTimePoint,nSrfs
52159599516SKenneth E. Jansen      real*8   pressHist(0:MAXSURF)
52259599516SKenneth E. Jansen      real*8   QHist(nTimePoint+1,nSrfs),betas(nTimePoint+2,nSrfs)
52359599516SKenneth E. Jansen      !don't need here betas(ntimePoint+2)
52459599516SKenneth E. Jansen      !but pb of array passing if cut at nTimePoint+1
52559599516SKenneth E. Jansen      pressHist=zero
52659599516SKenneth E. Jansen      do k=1,nSrfs
52759599516SKenneth E. Jansen        do j=1,nTimePoint+1
52859599516SKenneth E. Jansen            pressHist(k) = pressHist(k) + QHist(j,k)*betas(j,k)
52959599516SKenneth E. Jansen        enddo
53059599516SKenneth E. Jansen      enddo
53159599516SKenneth E. Jansen      return
53259599516SKenneth E. Jansen      end
533