1*59599516SKenneth E. Jansen subroutine genini (iBC, BC, y, ac, iper, ilwork, 2*59599516SKenneth E. Jansen & ifath, velbar, nsons,x, 3*59599516SKenneth E. Jansen & shp, shgl, shpb, shglb ) 4*59599516SKenneth E. Jansenc 5*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 6*59599516SKenneth E. Jansenc This routine reads the initial values in primitive form (density, 7*59599516SKenneth E. Jansenc velocity and temperature), satisfies the boundary conditions and 8*59599516SKenneth E. Jansenc converts them to Y-variables. 9*59599516SKenneth E. Jansenc 10*59599516SKenneth E. Jansenc input: 11*59599516SKenneth E. Jansenc iBC (nshg) : boundary condition code 12*59599516SKenneth E. Jansenc BC (nshg,ndofBC) : boundary condition constrain data 13*59599516SKenneth E. Jansenc x (numnp,nsd) : locations of nodes, numnp-> # of node 14*59599516SKenneth E. Jansenc nsd-> space dimension, 1=x, 2=y, 3=z 15*59599516SKenneth E. JansenC 16*59599516SKenneth E. Jansenc 17*59599516SKenneth E. Jansenc output: 18*59599516SKenneth E. Jansenc y (nshg,ndof) : initial values of Y variables 19*59599516SKenneth E. Jansenc 20*59599516SKenneth E. Jansenc 21*59599516SKenneth E. Jansenc Farzin Shakib, Winter 1986. 22*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991. (Fortran 90) 23*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 24*59599516SKenneth E. Jansenc 25*59599516SKenneth E. Jansen use specialBC ! gets itvn from here 26*59599516SKenneth E. Jansen use convolImpFlow ! brings in ntimeptpT and other variables 27*59599516SKenneth E. Jansen use convolRCRFlow ! brings RCR variables 28*59599516SKenneth E. Jansen include "common.h" 29*59599516SKenneth E. Jansen include "mpif.h" 30*59599516SKenneth E. Jansen include "auxmpi.h" 31*59599516SKenneth E. Jansenc 32*59599516SKenneth E. Jansen dimension iBC(nshg), iper(nshg), 33*59599516SKenneth E. Jansen & BC(nshg,ndofBC), y(nshg,ndof), 34*59599516SKenneth E. Jansen & ac(nshg,ndof), x(numnp,nsd), 35*59599516SKenneth E. Jansen & shp(MAXTOP,maxsh,MAXQPT), 36*59599516SKenneth E. Jansen & shgl(MAXTOP,nsd,maxsh,MAXQPT), 37*59599516SKenneth E. Jansen & shpb(MAXTOP,maxsh,MAXQPT), 38*59599516SKenneth E. Jansen & shglb(MAXTOP,nsd,maxsh,MAXQPT) 39*59599516SKenneth E. Jansen 40*59599516SKenneth E. Jansenc 41*59599516SKenneth E. Jansen dimension ilwork(nlwork) 42*59599516SKenneth E. Jansenc 43*59599516SKenneth E. Jansenc stuff for dynamic model s.w.avg and wall model 44*59599516SKenneth E. Jansenc 45*59599516SKenneth E. Jansen dimension ifath(numnp), velbar(nfath,nflow), 46*59599516SKenneth E. Jansen & nsons(nfath) 47*59599516SKenneth E. Jansen 48*59599516SKenneth E. Jansen character*20 fname1 49*59599516SKenneth E. Jansen character*10 cname2 50*59599516SKenneth E. Jansen character*5 cname 51*59599516SKenneth E. Jansenc 52*59599516SKenneth E. Jansenc.... --------------------------> Restart <--------------------------- 53*59599516SKenneth E. Jansenc 54*59599516SKenneth E. Jansenc.... read q from [RESTAR.INP], reset LSTEP 55*59599516SKenneth E. Jansenc 56*59599516SKenneth E. Jansen call restar ('in ', y, ac) 57*59599516SKenneth E. Jansenc 58*59599516SKenneth E. Jansen if((itwmod.gt.0) 59*59599516SKenneth E. Jansen & .or. (nsonmax.eq.1 .and. iLES.gt.0) ) then 60*59599516SKenneth E. Jansen call rwvelb('in ',velbar,ifail) 61*59599516SKenneth E. Jansenc 62*59599516SKenneth E. Jansenc if the read failed calculate velbar 63*59599516SKenneth E. Jansenc 64*59599516SKenneth E. Jansen if(ifail.eq.1) then 65*59599516SKenneth E. Jansen call getvel (y, ilwork, iBC, 66*59599516SKenneth E. Jansen & nsons, ifath, velbar) 67*59599516SKenneth E. Jansen endif 68*59599516SKenneth E. Jansen 69*59599516SKenneth E. Jansen endif ! for the itwmod or irscale 70*59599516SKenneth E. Jansenc 71*59599516SKenneth E. Jansenc.... time varying boundary conditions as set from file bct.dat and impt.dat 72*59599516SKenneth E. Jansenc (see function for format in file bctint.f) 73*59599516SKenneth E. Jansenc 74*59599516SKenneth E. Jansen if (itvn .gt. 0 ) then !for inlet velocities 75*59599516SKenneth E. Jansen call initBCt( x, iBC, BC) 76*59599516SKenneth E. Jansen call BCint(lstep*Delt(1),shp,shgl,shpb,shglb,x, BC, iBC) 77*59599516SKenneth E. Jansen endif 78*59599516SKenneth E. Jansen if (impfile .gt. 0 ) then !for impedance BC 79*59599516SKenneth E. Jansen if(myrank.eq.master) then 80*59599516SKenneth E. Jansen write(*,*) 'reading Qhistor.dat' 81*59599516SKenneth E. Jansen endif 82*59599516SKenneth E. Jansen open(unit=816, file='Qhistor.dat',status='old') 83*59599516SKenneth E. Jansen read (816,*) ntimeptpT 84*59599516SKenneth E. Jansen allocate (QHistImp(ntimeptpT+1,numImpSrfs)) 85*59599516SKenneth E. Jansen allocate (QHistTry(ntimeptpT,numImpSrfs)) 86*59599516SKenneth E. Jansen allocate (QHistTryF(ntimeptpT,numImpSrfs)) 87*59599516SKenneth E. Jansen do j=1,ntimeptpT+1 88*59599516SKenneth E. Jansen read(816,*) (QHistImp(j,n),n=1,numImpSrfs) !read flow history 89*59599516SKenneth E. Jansen enddo 90*59599516SKenneth E. Jansen close(816) 91*59599516SKenneth E. Jansen call initImpt() !read impedance data and initialize begin/end values 92*59599516SKenneth E. Jansen do i=2,ntimeptpT 93*59599516SKenneth E. Jansen call Impint((ntimeptpT-i+1)*Delt(1),i) !return Imp values in reverse order ZN->Z0 94*59599516SKenneth E. Jansen enddo 95*59599516SKenneth E. Jansen 96*59599516SKenneth E. Jansen allocate (poldImp(0:MAXSURF)) !for pressure part that depends on the history only 97*59599516SKenneth E. Jansen endif 98*59599516SKenneth E. Jansen if (ircrfile .gt. 0 ) then !for RCR BC 99*59599516SKenneth E. Jansen call initRCRt() !read RCR data 100*59599516SKenneth E. Jansen dtRCR(:) = Delt(1)/(ValueListRCR(2,:)*ValueListRCR(3,:)) 101*59599516SKenneth E. Jansen allocate (RCRConvCoef(nstep(1)+lstep+2,numRCRSrfs)) !for convolution coeff 102*59599516SKenneth E. Jansen !last one needed only to have array of same size as imp BC 103*59599516SKenneth E. Jansen allocate (QHistRCR(nstep(1)+lstep+1,numRCRSrfs)) !for flow history 104*59599516SKenneth E. Jansen QHistRCR = zero 105*59599516SKenneth E. Jansen if(lstep.ne.0) then 106*59599516SKenneth E. Jansen fname1 = 'Qhistor' 107*59599516SKenneth E. Jansen fname1 = trim(fname1)//trim(cname2(lstep))//'.dat' 108*59599516SKenneth E. Jansen fname1 = trim(fname1) 109*59599516SKenneth E. Jansen if(myrank.eq.master) then 110*59599516SKenneth E. Jansen write(*,*) 'reading ', fname1 111*59599516SKenneth E. Jansen endif 112*59599516SKenneth E. Jansen open(unit=816, file=fname1, status='old') 113*59599516SKenneth E. Jansen read(816,*) k 114*59599516SKenneth E. Jansen do j=1,lstep+1 115*59599516SKenneth E. Jansen read(816,*) (QHistRCR(j,n),n=1,numRCRSrfs) !read flow history 116*59599516SKenneth E. Jansen enddo 117*59599516SKenneth E. Jansen close(816) 118*59599516SKenneth E. Jansen endif 119*59599516SKenneth E. Jansen allocate (poldRCR(0:MAXSURF)) !for pressure part that depends on the history only 120*59599516SKenneth E. Jansen allocate (HopRCR(0:MAXSURF)) !for H operator contribution 121*59599516SKenneth E. Jansen endif 122*59599516SKenneth E. Jansen 123*59599516SKenneth E. Jansen 124*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc 125*59599516SKenneth E. Jansenc this subroutine initBCprofileScale is called to generate the correct 126*59599516SKenneth E. Jansenc BC array, including the siuation of Take BC from IC for inlet. 127*59599516SKenneth E. Jansenc so this subroutine must be called before BCs are applied to ICs 128*59599516SKenneth E. Jansenc as those following lines do 129*59599516SKenneth E. Jansenc call INIBCprofile(BC,y,x) 130*59599516SKenneth E. Jansenc call MPI_BARRIER(MPI_COMM_WORLD,ierr) 131*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc 132*59599516SKenneth E. Jansenc 133*59599516SKenneth E. Jansenc.... satisfy the boundary conditions 134*59599516SKenneth E. Jansenc 135*59599516SKenneth E. Jansen call itrBC (y, ac, iBC, BC, iper, ilwork) 136*59599516SKenneth E. Jansen 137*59599516SKenneth E. Jansen itempr=mod(impl(1),2) ! tempr solve if impl odd 138*59599516SKenneth E. Jansen if(itempr.eq.1) then 139*59599516SKenneth E. Jansen isclr=0 140*59599516SKenneth E. Jansen call itrBCSclr (y, ac, iBC, BC, iper, ilwork) 141*59599516SKenneth E. Jansen endif 142*59599516SKenneth E. Jansen do isclr=1,nsclr 143*59599516SKenneth E. Jansen call itrBCSclr (y, ac, iBC, BC, iper, ilwork) 144*59599516SKenneth E. Jansen enddo 145*59599516SKenneth E. Jansen 146*59599516SKenneth E. Jansen if((irscale.ge.0) .and. (myrank.eq.master)) then 147*59599516SKenneth E. Jansen call genscale(y, x, iBC) 148*59599516SKenneth E. Jansen endif 149*59599516SKenneth E. Jansenc 150*59599516SKenneth E. Jansenc.... ---------------------------> Echo <---------------------------- 151*59599516SKenneth E. Jansenc 152*59599516SKenneth E. Jansenc.... echo the initial data 153*59599516SKenneth E. Jansenc 154*59599516SKenneth E. Jansen if ((necho .lt. 0).and.(myrank.eq.master)) then 155*59599516SKenneth E. Jansen do n = 1, nshg 156*59599516SKenneth E. Jansen if (mod(n,50) .eq. 1) write(iecho,1000) ititle,(i,i=1,ndof) 157*59599516SKenneth E. Jansen write (iecho,1100) n, (y(n,i),i=1,ndof) 158*59599516SKenneth E. Jansen enddo 159*59599516SKenneth E. Jansen endif 160*59599516SKenneth E. Jansenc 161*59599516SKenneth E. Jansenc.... return 162*59599516SKenneth E. Jansenc 163*59599516SKenneth E. Jansen return 164*59599516SKenneth E. Jansenc 165*59599516SKenneth E. Jansen1000 format(a80,//, 166*59599516SKenneth E. Jansen & ' I n i t i a l V a l u e s ',//, 167*59599516SKenneth E. Jansen & ' Node ',/, 168*59599516SKenneth E. Jansen & ' Number ',6x,6('dof',i1,:,10x)) 169*59599516SKenneth E. Jansen1100 format(1p,2x,i5,5x,5(e12.5,2x)) 170*59599516SKenneth E. Jansenc 171*59599516SKenneth E. Jansen end 172