xref: /phasta/phSolver/compressible/MachControl.f90 (revision 1016729149754f57cd03fe576ba6fd0f1723ab31)
1*10167291SKenneth E. Jansen      module MachControl
2*10167291SKenneth E. Jansen        logical exMC            !Mach Control enable / Does the MC file
3*10167291SKenneth E. Jansen                                !exist
4*10167291SKenneth E. Jansen
5*10167291SKenneth E. Jansen        integer :: MC_iPP       !Probe Point index
6*10167291SKenneth E. Jansen        integer :: MC_surfID    !Surface ID to identify where to apply
7*10167291SKenneth E. Jansen                                ! the boundary condition
8*10167291SKenneth E. Jansen
9*10167291SKenneth E. Jansen        real*8 :: MC_dt         !simulaiton time step
10*10167291SKenneth E. Jansen        real*8 :: MC_MtTarget   !Target throat Mach number
11*10167291SKenneth E. Jansen        real*8 :: MC_dMtCO      !M_throat Cut Off (frequency) (for
12*10167291SKenneth E. Jansen                                ! filtering the state)
13*10167291SKenneth E. Jansen        real*8 :: MC_ulim(2)    !u_inlet velocity limit
14*10167291SKenneth E. Jansen                                !u_inlet(1)     lower limit
15*10167291SKenneth E. Jansen                                !u_inlet(2)     upper limit
16*10167291SKenneth E. Jansen        !States
17*10167291SKenneth E. Jansen        real*8 :: MC_Mt         !actual Mach number
18*10167291SKenneth E. Jansen        real*8 :: MC_dMtdt      !time derivative of MC_dMt
19*10167291SKenneth E. Jansen        real*8 :: MC_dMt        !filtered state deviation from target
20*10167291SKenneth E. Jansen        real*8 :: MC_IdMt       !integreal of MC_dMt
21*10167291SKenneth E. Jansen
22*10167291SKenneth E. Jansen        !Gains
23*10167291SKenneth E. Jansen        real*8 :: MC_Kd         !derivative gain
24*10167291SKenneth E. Jansen        real*8 :: MC_Kp         !proportional gain
25*10167291SKenneth E. Jansen        real*8 :: MC_KI         !integral gain
26*10167291SKenneth E. Jansen
27*10167291SKenneth E. Jansen        !Variables to apply the boundary condition
28*10167291SKenneth E. Jansen        integer, allocatable, dimension(:) :: MC_imapInlet(:)
29*10167291SKenneth E. Jansen        integer :: MC_icountInlet   !number of inlet boundary points
30*10167291SKenneth E. Jansen        real*8, allocatable, dimension(:) :: MC_vscale(:)   !scaling factor for BL
31*10167291SKenneth E. Jansen        real*8 :: MC_BLdelta        !boundary layer thickness
32*10167291SKenneth E. Jansen        real*8 :: MC_uInlet     !x-velocity at the inlet
33*10167291SKenneth E. Jansen
34*10167291SKenneth E. Jansen      end module
35*10167291SKenneth E. Jansen
36*10167291SKenneth E. Jansen      subroutine MC_init(dt, tstep, BC)
37*10167291SKenneth E. Jansen        use MachControl
38*10167291SKenneth E. Jansen        use turbSA
39*10167291SKenneth E. Jansen        include "common.h"
40*10167291SKenneth E. Jansen        include "mpif.h"
41*10167291SKenneth E. Jansen
42*10167291SKenneth E. Jansen        dimension :: BC(nshg, ndofBC)
43*10167291SKenneth E. Jansen        logical :: MC_existDir
44*10167291SKenneth E. Jansen        real*8 :: dt
45*10167291SKenneth E. Jansen        real*8 :: vBC, vBCg
46*10167291SKenneth E. Jansen        integer :: tstep
47*10167291SKenneth E. Jansen        vBCg = 0
48*10167291SKenneth E. Jansen
49*10167291SKenneth E. Jansen        if(myrank == master) then
50*10167291SKenneth E. Jansen          inquire(file="./MachControl/.", exist=MC_existDir)
51*10167291SKenneth E. Jansen          if(.not. MC_existDir) then
52*10167291SKenneth E. Jansen            call system("mkdir ./MachControl")  !Doesn't seem to work on BGQ.
53*10167291SKenneth E. Jansen          endif
54*10167291SKenneth E. Jansen        endif
55*10167291SKenneth E. Jansen
56*10167291SKenneth E. Jansen        call MC_readState(tstep)
57*10167291SKenneth E. Jansen
58*10167291SKenneth E. Jansen        if(exMC) then
59*10167291SKenneth E. Jansen          !find the node associated with isetBlowerID_Duct and write the
60*10167291SKenneth E. Jansen          !scale vector to linear ramp velocity near the walls.
61*10167291SKenneth E. Jansen          allocate(MC_imapInlet(nshg))
62*10167291SKenneth E. Jansen          allocate(MC_vscale(nshg))
63*10167291SKenneth E. Jansen          call sfID2np(MC_surfID, MC_icountInlet, MC_imapInlet)
64*10167291SKenneth E. Jansen
65*10167291SKenneth E. Jansen          if(MC_BLdelta > 0) then
66*10167291SKenneth E. Jansen            do i = 1, MC_icountInlet
67*10167291SKenneth E. Jansen              imapped = MC_imapInlet(i)
68*10167291SKenneth E. Jansen
69*10167291SKenneth E. Jansen              MC_vscale(i) = min(one, d2wall(imapped)/MC_BLdelta)
70*10167291SKenneth E. Jansen              !MC_vscale(i) = d2wall(imapped)/MC_BLdelta
71*10167291SKenneth E. Jansen              !if(MC_vscale(i) > 1) MC_vscale(i) = 1
72*10167291SKenneth E. Jansen
73*10167291SKenneth E. Jansen              !Also look for the inlet velocity. This will get
74*10167291SKenneth E. Jansen              !overwritten in the next time step, but it's nice to know
75*10167291SKenneth E. Jansen              !what it started with.
76*10167291SKenneth E. Jansen              vBCg = max(vBCg, BC(imapped, 3))
77*10167291SKenneth E. Jansen            enddo
78*10167291SKenneth E. Jansen          endif
79*10167291SKenneth E. Jansen
80*10167291SKenneth E. Jansen          !Look for the maximum global inlet velocity among all parts.
81*10167291SKenneth E. Jansen          if(numpe > 1) then
82*10167291SKenneth E. Jansen            call MPI_ALLREDUCE(vBCg, MC_uInlet, 1,
83*10167291SKenneth E. Jansen     &           MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr)
84*10167291SKenneth E. Jansen          endif
85*10167291SKenneth E. Jansen
86*10167291SKenneth E. Jansen          MC_dt = dt
87*10167291SKenneth E. Jansen        endif
88*10167291SKenneth E. Jansen      end subroutine
89*10167291SKenneth E. Jansen
90*10167291SKenneth E. Jansen
91*10167291SKenneth E. Jansen      subroutine MC_updateState()
92*10167291SKenneth E. Jansen        !Updates the derivative, proportional and integral values used
93*10167291SKenneth E. Jansen        !when calculating MC_uInlet.
94*10167291SKenneth E. Jansen        !Note: This function uses the most recent values in the TimeData
95*10167291SKenneth E. Jansen        !buffer vartsbuff. TD_bufferData should be called prior to this
96*10167291SKenneth E. Jansen        !subroutine.
97*10167291SKenneth E. Jansen
98*10167291SKenneth E. Jansen        use MachControl
99*10167291SKenneth E. Jansen        use TimeData
100*10167291SKenneth E. Jansen        include "common.h"
101*10167291SKenneth E. Jansen        include "mpif.h"
102*10167291SKenneth E. Jansen
103*10167291SKenneth E. Jansen        real*8 :: u2, T
104*10167291SKenneth E. Jansen        real*8 :: Mt, dMt, MC_alpha
105*10167291SKenneth E. Jansen        logical :: bangbang
106*10167291SKenneth E. Jansen
107*10167291SKenneth E. Jansen        if(myrank == master) then !MC_* are not set on other processors
108*10167291SKenneth E. Jansen          T  = vartsbuff(MC_iPP, 5, ivartsBuff)
109*10167291SKenneth E. Jansen          u2 = vartsbuff(MC_iPP, 2, ivartsBuff)**2 +
110*10167291SKenneth E. Jansen     &         vartsbuff(MC_iPP, 3, ivartsBuff)**2 +
111*10167291SKenneth E. Jansen     &         vartsbuff(MC_iPP, 4, ivartsBuff)**2
112*10167291SKenneth E. Jansen
113*10167291SKenneth E. Jansen          MC_Mt = sqrt(u2/(1.4*Rgas*T))
114*10167291SKenneth E. Jansen          dMt = MC_Mt - MC_MtTarget
115*10167291SKenneth E. Jansen
116*10167291SKenneth E. Jansen          !Note: MC_uInlet should always equal the previous time step;
117*10167291SKenneth E. Jansen          !on the first time step of a run, MC_uInlet is computed based
118*10167291SKenneth E. Jansen          !off of the maximum inlet velocity in the solution.
119*10167291SKenneth E. Jansen          bangbang = MC_uInlet <= MC_uLim(1) .or.
120*10167291SKenneth E. Jansen     &               MC_uInlet >= MC_uLim(2)
121*10167291SKenneth E. Jansen
122*10167291SKenneth E. Jansen          !Low pass filter dMt
123*10167291SKenneth E. Jansen          MC_alpha = MC_dt/(MC_dMtCO + MC_dt)
124*10167291SKenneth E. Jansen          dMt = MC_alpha*dMt + (1 - MC_alpha)*MC_dMt
125*10167291SKenneth E. Jansen
126*10167291SKenneth E. Jansen          !Calculate derivatives and integrals
127*10167291SKenneth E. Jansen          if(.not. bangbang) then
128*10167291SKenneth E. Jansen            !Only integrate when not hitting the limits. This introduces
129*10167291SKenneth E. Jansen            !a nonlinearity in the control scheme which minimizes the
130*10167291SKenneth E. Jansen            !overshoot when starting from zero initial conditions
131*10167291SKenneth E. Jansen            MC_IdMt  = (dMt + MC_dMt)*MC_dt/2 + MC_IdMt !trapezoidal rule
132*10167291SKenneth E. Jansen          endif
133*10167291SKenneth E. Jansen          MC_dMtdt = (dMt - MC_dMt)/MC_dt !first order BE FD derivative
134*10167291SKenneth E. Jansen          MC_dMt = dMt !update the actual state
135*10167291SKenneth E. Jansen        endif
136*10167291SKenneth E. Jansen      end subroutine
137*10167291SKenneth E. Jansen
138*10167291SKenneth E. Jansen
139*10167291SKenneth E. Jansen      subroutine MC_applyBC(BC)
140*10167291SKenneth E. Jansen        use MachControl
141*10167291SKenneth E. Jansen        include "common.h"
142*10167291SKenneth E. Jansen        include "mpif.h"
143*10167291SKenneth E. Jansen
144*10167291SKenneth E. Jansen        dimension :: BC(nshg, ndofBC)
145*10167291SKenneth E. Jansen        real*8 :: vInlet
146*10167291SKenneth E. Jansen        integer :: imapped
147*10167291SKenneth E. Jansen        real*8 :: tmp
148*10167291SKenneth E. Jansen
149*10167291SKenneth E. Jansen        if(myrank .eq. master) then
150*10167291SKenneth E. Jansen          MC_uInlet = MC_Kd*MC_dMtdt +
151*10167291SKenneth E. Jansen     &                MC_Kp*MC_dMt +
152*10167291SKenneth E. Jansen     &                MC_KI*MC_IdMt
153*10167291SKenneth E. Jansen
154*10167291SKenneth E. Jansen          !clip to [MC_uLim(1),MC_uLim(2)]
155*10167291SKenneth E. Jansen          if(MC_uLim(1) < MC_ulim(2)) then
156*10167291SKenneth E. Jansen            MC_uInlet = min(max(MC_uInlet, MC_uLim(1)), MC_uLim(2))
157*10167291SKenneth E. Jansen          endif
158*10167291SKenneth E. Jansen        endif
159*10167291SKenneth E. Jansen        if(numpe .gt. 1)
160*10167291SKenneth E. Jansen     &    call MPI_Bcast(MC_uInlet, 1, MPI_DOUBLE_PRECISION,
161*10167291SKenneth E. Jansen     &                         master, MPI_COMM_WORLD, ierr)
162*10167291SKenneth E. Jansen
163*10167291SKenneth E. Jansen        if(MC_BLdelta .gt. 0) then !Apply a velocity profile with a BL
164*10167291SKenneth E. Jansen          do i = 1, MC_icountInlet
165*10167291SKenneth E. Jansen            imapped = MC_imapInlet(i)
166*10167291SKenneth E. Jansen            BC(imapped,3) = MC_uInlet*MC_vscale(i)
167*10167291SKenneth E. Jansen          enddo
168*10167291SKenneth E. Jansen        else   !Apply a uniform inlet velocity.
169*10167291SKenneth E. Jansen          do i = 1, MC_icountInlet
170*10167291SKenneth E. Jansen            BC(MC_imapInlet(i), 3) = MC_uInlet
171*10167291SKenneth E. Jansen          enddo
172*10167291SKenneth E. Jansen        endif
173*10167291SKenneth E. Jansen
174*10167291SKenneth E. Jansen      end subroutine
175*10167291SKenneth E. Jansen
176*10167291SKenneth E. Jansen
177*10167291SKenneth E. Jansen      subroutine MC_printState()
178*10167291SKenneth E. Jansen      !Prints the Mach Control states into the output
179*10167291SKenneth E. Jansen
180*10167291SKenneth E. Jansen        use MachControl
181*10167291SKenneth E. Jansen        include "common.h"
182*10167291SKenneth E. Jansen
183*10167291SKenneth E. Jansen        if(myrank == master) then
184*10167291SKenneth E. Jansen          write(*, "('M_throat: ',F8.6,'  dM_throat: ',F8.6)")
185*10167291SKenneth E. Jansen     &             MC_Mt, MC_Mt - MC_MtTarget
186*10167291SKenneth E. Jansen          write(*,"('dM/dt = ',E14.6,'  M  = ',E14.6,'  I(M) =',E14.6)")
187*10167291SKenneth E. Jansen     &             MC_dMtdt, MC_dMt, MC_IdMt
188*10167291SKenneth E. Jansen          write(*,"('Kd    = ',E14.6,'  Kp = ',E14.6,'  KI   =',E14.6)")
189*10167291SKenneth E. Jansen     &            MC_Kd, MC_Kp, MC_KI
190*10167291SKenneth E. Jansen          write(*, "('u_inlet: ',F9.6)") MC_uInlet
191*10167291SKenneth E. Jansen        endif
192*10167291SKenneth E. Jansen      end subroutine
193*10167291SKenneth E. Jansen
194*10167291SKenneth E. Jansen      subroutine MC_readState(tstep)
195*10167291SKenneth E. Jansen        !Reads the intput / state file for MachControl.
196*10167291SKenneth E. Jansen        !Written by: Nicholas Mati      2014-04-19
197*10167291SKenneth E. Jansen        !Revision History:
198*10167291SKenneth E. Jansen        ! - 2014-04-19: created
199*10167291SKenneth E. Jansen        !
200*10167291SKenneth E. Jansen        !The file is expected to be of the form:
201*10167291SKenneth E. Jansen        ! MC_iPP   MC_surfID  MC_MtTarget
202*10167291SKenneth E. Jansen        ! MC_dMtdt MC_dMt MC_IdMt
203*10167291SKenneth E. Jansen        ! MC_Kd    MC_Kp  MC_KI
204*10167291SKenneth E. Jansen        ! MC_dMtCO MC_BLdelta
205*10167291SKenneth E. Jansen
206*10167291SKenneth E. Jansen        use MachControl
207*10167291SKenneth E. Jansen        include "common.h"
208*10167291SKenneth E. Jansen        include "mpif.h"
209*10167291SKenneth E. Jansen
210*10167291SKenneth E. Jansen        integer :: tstep
211*10167291SKenneth E. Jansen        character(len=60) :: fname
212*10167291SKenneth E. Jansen
213*10167291SKenneth E. Jansen        if(myrank == master) then
214*10167291SKenneth E. Jansen          write(fname, "('./MachControl/MachControl.',I0,'.dat')") tstep
215*10167291SKenneth E. Jansen          inquire(file=fname, exist=exMC)
216*10167291SKenneth E. Jansen
217*10167291SKenneth E. Jansen          if(exMC) then
218*10167291SKenneth E. Jansen            open(unit=1002, file=fname, status='old')
219*10167291SKenneth E. Jansen
220*10167291SKenneth E. Jansen            read(1002, *)  MC_iPP,   MC_surfID
221*10167291SKenneth E. Jansen            read(1002, *)  MC_Mt,    MC_MtTarget
222*10167291SKenneth E. Jansen            read(1002, *)  MC_dMtdt, MC_dMt,    MC_IdMt
223*10167291SKenneth E. Jansen            read(1002, *)  MC_Kd,    MC_Kp,     MC_KI
224*10167291SKenneth E. Jansen            read(1002, *)  MC_dMtCO, MC_BLdelta
225*10167291SKenneth E. Jansen            read(1002, *)  MC_uLim(1), MC_uLim(2)
226*10167291SKenneth E. Jansen            close(1002)
227*10167291SKenneth E. Jansen
228*10167291SKenneth E. Jansen            !If d2wall is not available, then set the inlet BL thickness
229*10167291SKenneth E. Jansen            !to zero to disable wall scalling.
230*10167291SKenneth E. Jansen            if(iRANS >= 0) MC_BLdelta = 0
231*10167291SKenneth E. Jansen          endif
232*10167291SKenneth E. Jansen        endif
233*10167291SKenneth E. Jansen
234*10167291SKenneth E. Jansen        !Only a few values are needed across other processors.
235*10167291SKenneth E. Jansen         if(numpe > 1) then
236*10167291SKenneth E. Jansen          call MPI_BCAST(MC_BLdelta, 1, MPI_DOUBLE_PRECISION,
237*10167291SKenneth E. Jansen     &                          master, MPI_COMM_WORLD, ierr)
238*10167291SKenneth E. Jansen          call MPI_BCAST(MC_surfID,  1, MPI_INTEGER,
239*10167291SKenneth E. Jansen     &                          master, MPI_COMM_WORLD, ierr)
240*10167291SKenneth E. Jansen          call MPI_BCAST(exMC,       1, MPI_INTEGER,
241*10167291SKenneth E. Jansen     &                          master, MPI_COMM_WORLD, ierr)
242*10167291SKenneth E. Jansen        endif
243*10167291SKenneth E. Jansen
244*10167291SKenneth E. Jansen      end subroutine
245*10167291SKenneth E. Jansen
246*10167291SKenneth E. Jansen      subroutine MC_writeState(tstep)
247*10167291SKenneth E. Jansen        !Writes a file with the name MachControl.[ts].dat that can later
248*10167291SKenneth E. Jansen        !be read in and restarted from.
249*10167291SKenneth E. Jansen
250*10167291SKenneth E. Jansen        use MachControl
251*10167291SKenneth E. Jansen        include "common.h"
252*10167291SKenneth E. Jansen
253*10167291SKenneth E. Jansen        integer tstep
254*10167291SKenneth E. Jansen        character(len=60) :: fname
255*10167291SKenneth E. Jansen
256*10167291SKenneth E. Jansen        if(myrank == master) then
257*10167291SKenneth E. Jansen          write(fname, "('./MachControl/MachControl.',I0,'.dat')") tstep
258*10167291SKenneth E. Jansen          open(unit=1002, file=fname, status='unknown')
259*10167291SKenneth E. Jansen
260*10167291SKenneth E. Jansen          write(1002, "(2(i9))")         MC_iPP,   MC_surfID
261*10167291SKenneth E. Jansen          write(1002, "(2(F18.14))")     MC_Mt,    MC_MtTarget
262*10167291SKenneth E. Jansen          write(1002, "(3(F18.14))")     MC_dMtdt, MC_dMt,    MC_IdMt
263*10167291SKenneth E. Jansen          write(1002, "(3(F18.10))")     MC_Kd,    MC_Kp,     MC_KI
264*10167291SKenneth E. Jansen          write(1002, "(2(F18.14))")     MC_dMtCO, MC_BLdelta
265*10167291SKenneth E. Jansen          write(1002, "(2(F18.14))")     MC_uLim(1), MC_uLim(2)
266*10167291SKenneth E. Jansen          close(1002)
267*10167291SKenneth E. Jansen        endif
268*10167291SKenneth E. Jansen      end subroutine
269*10167291SKenneth E. Jansen
270*10167291SKenneth E. Jansen      subroutine MC_finalize()
271*10167291SKenneth E. Jansen
272*10167291SKenneth E. Jansen        use MachControl
273*10167291SKenneth E. Jansen
274*10167291SKenneth E. Jansen        deallocate(MC_imapInlet)
275*10167291SKenneth E. Jansen        deallocate(MC_vscale)
276*10167291SKenneth E. Jansen      end subroutine
277*10167291SKenneth E. Jansen
278*10167291SKenneth E. Jansen
279