Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLO_ModBoundaryConditions.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 ! ******************************************************************************
24 !
25 ! Purpose: Suite of boundary condition routines.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLO_ModBoundaryConditions.F90,v 1.30 2010/02/18 21:47:39 juzhang Exp $
34 !
35 ! Copyright: (c) 2004 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY: t_global
42  USE modparameters
43  USE moddatatypes
44  USE moderror
45  USE moddatastruct, ONLY: t_region
46  USE modbndpatch, ONLY : t_patch
47  USE modgrid, ONLY: t_grid
48  USE modmpi
49 
50  IMPLICIT NONE
51 
52  PRIVATE
53  PUBLIC :: rflo_bcondfarfield, &
65 
66 ! private : RFLO_BcondInflowVTPerf
67 ! RFLO_BcondInflowVPPerf
68 ! RFLO_BcondInflowFluc
69 ! RFLO_BcondInjectionAPN
70 
71 ! ******************************************************************************
72 ! Declarations and definitions
73 ! ******************************************************************************
74 
75  CHARACTER(CHRLEN) :: RCSIdentString = &
76  '$RCSfile: RFLO_ModBoundaryConditions.F90,v $ $Revision: 1.30 $'
77 
78 ! ******************************************************************************
79 ! Routines
80 ! ******************************************************************************
81 
82  CONTAINS
83 
84 !******************************************************************************
85 !
86 ! Purpose: update values in dummy cells at far field boundary patch.
87 !
88 ! Description: none.
89 !
90 ! Input: region = region dimensions, user input
91 ! patch = current patch.
92 !
93 ! Output: region%levels%mixt = flow variables in dummy cells.
94 !
95 ! Notes: there is no vortex correction implemented yet.
96 !
97 !******************************************************************************
98 
99 SUBROUTINE rflo_bcondfarfield( region,patch )
100 
104  IMPLICIT NONE
105 
106 #include "Indexing.h"
107 
108 ! ... parameters
109  TYPE(t_region) :: region
110  TYPE(t_patch) :: patch
111 
112 ! ... loop variables
113  INTEGER :: idum, i, j, k, n1, n2
114 
115 ! ... local variables
116  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
117  INTEGER :: ilev, lbound, icoff, ijcoff, inoff, ijnoff, noff, distrib
118  INTEGER :: indcp, indmol, gasmodel, ijkc, ijkc1, ijkd, ijkn, i2d
119  INTEGER :: inode, jnode, knode
120 
121  REAL(RFREAL) :: sgn, ds, sxn, syn, szn
122  REAL(RFREAL) :: mach, attack, slip, press, temp, pb
123  REAL(RFREAL), POINTER :: cv(:,:), dv(:,:), gv(:,:), vals(:,:), sface(:,:)
124 
125 !******************************************************************************
126 
127  CALL registerfunction( region%global,'RFLO_BcondFarfield',&
128  'RFLO_ModBoundaryConditions.F90' )
129 
130 ! get dimensions and pointers
131 
132  ilev = region%currLevel
133  lbound = patch%lbound
134 
135  CALL rflo_getpatchindices( region,patch,ilev, &
136  ibeg,iend,jbeg,jend,kbeg,kend )
138  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
139  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
140 
141  noff = abs(patch%l1end-patch%l1beg) + 1
142  distrib = patch%mixt%distrib
143  indcp = region%levels(ilev)%mixt%indCp
144  indmol = region%levels(ilev)%mixt%indMol
145  gasmodel = region%mixtInput%gasModel
146 
147  cv => region%levels(ilev)%mixt%cv
148  dv => region%levels(ilev)%mixt%dv
149  gv => region%levels(ilev)%mixt%gv
150  vals => patch%mixt%vals
151 
152 ! to take the right face vector and make it point outwards
153 
154  sgn = +1._rfreal
155  inode = 0
156  jnode = 0
157  knode = 0
158  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
159  sgn = -1._rfreal
160  inode = -idir
161  jnode = -jdir
162  knode = -kdir
163  ENDIF
164 
165 ! get the appropriate face vector
166 
167  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
168  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
169  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
170 
171 ! loop over all cells of a patch
172 
173  DO idum=1,region%nDumCells
174  DO k=kbeg,kend
175  DO j=jbeg,jend
176  DO i=ibeg,iend
177  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
178 
179  IF (idum == 1) THEN
180  ijkc = indijk(i,j,k,icoff,ijcoff)
181  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
182  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
183  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
184  sface(zcoord,ijkn)*sface(zcoord,ijkn))
185  sxn = sgn*sface(xcoord,ijkn)/ds
186  syn = sgn*sface(ycoord,ijkn)/ds
187  szn = sgn*sface(zcoord,ijkn)/ds
188 
189  IF (lbound==1 .OR. lbound==2) THEN
190  n1 = j - jbeg
191  n2 = k - kbeg
192  ELSE IF (lbound==3 .OR. lbound==4) THEN
193  n1 = k - kbeg
194  n2 = i - ibeg
195  ELSE IF (lbound==5 .OR. lbound==6) THEN
196  n1 = i - ibeg
197  n2 = j - jbeg
198  ENDIF
199  i2d = distrib * indij(n1,n2,noff)
200 
201  IF (gasmodel == gas_model_tcperf) THEN
202  mach = vals(bcdat_farf_mach ,i2d)
203  attack = vals(bcdat_farf_attack,i2d)
204  slip = vals(bcdat_farf_slip ,i2d)
205  press = vals(bcdat_farf_press ,i2d)
206  temp = vals(bcdat_farf_temp ,i2d)
207  CALL bcondfarfieldperf( mach,attack,slip,press,temp,sxn,syn,szn, &
208  gv(gv_mixt_cp ,ijkc*indcp),gv(gv_mixt_mol ,ijkc*indmol), &
209  cv(cv_mixt_dens,ijkc) ,cv(cv_mixt_xmom,ijkc) , &
210  cv(cv_mixt_ymom,ijkc) ,cv(cv_mixt_zmom,ijkc) , &
211  cv(cv_mixt_ener,ijkc) ,dv(dv_mixt_pres,ijkc) , &
212  cv(cv_mixt_dens,ijkd) ,cv(cv_mixt_xmom,ijkd) , &
213  cv(cv_mixt_ymom,ijkd) ,cv(cv_mixt_zmom,ijkd) , &
214  cv(cv_mixt_ener,ijkd) ,pb )
215  ELSE
216  CALL errorstop( region%global,err_unknown_bc,&
217  __line__ )
218  ENDIF
219 
220  ELSE ! idum > 1
221  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
222  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
223 
224  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
225  cv(cv_mixt_dens,ijkc1)
226  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
227  cv(cv_mixt_xmom,ijkc1)
228  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
229  cv(cv_mixt_ymom,ijkc1)
230  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
231  cv(cv_mixt_zmom,ijkc1)
232  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
233  cv(cv_mixt_ener,ijkc1)
234  ENDIF
235 
236  IF (gasmodel == gas_model_tcperf) THEN
237  CALL mixtureproperties( region,ijkd,ijkd,.false. )
238  ELSE
239  CALL mixtureproperties( region,ijkd,ijkd,.true. )
240  ENDIF
241 
242  ENDDO ! i
243  ENDDO ! j
244  ENDDO ! k
245  ENDDO ! idum
246 
247 ! finalize
248 
249  CALL deregisterfunction( region%global )
250 
251 END SUBROUTINE rflo_bcondfarfield
252 
253 !******************************************************************************
254 !
255 ! Purpose: update values in dummy cells at inflow boundary patch.
256 !
257 ! Description: none.
258 !
259 ! Input: region = region dimensions, user input
260 ! patch = current patch.
261 !
262 ! Output: region%levels%mixt = flow variables in dummy cells.
263 !
264 ! Notes: none.
265 !
266 !******************************************************************************
267 
268 SUBROUTINE rflo_bcondinflow( region,patch )
269 
273  IMPLICIT NONE
274 
275 #include "Indexing.h"
276 
277 ! ... parameters
278  TYPE(t_region) :: region
279  TYPE(t_patch) :: patch
280 
281 ! ... loop variables
282  INTEGER :: idum, i, j, k, n1, n2
283 
284 ! ... local variables
285  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
286  INTEGER :: ilev, lbound, icoff, ijcoff, inoff, ijnoff, noff, bctype, &
287  distrib, bcoptfixed, bcopttype, bcoptmodel
288  INTEGER :: indcp, indmol, gasmodel, ijkc, ijkc1, ijkd, ijkn, i2d
289  INTEGER :: inode, jnode, knode
290  INTEGER :: ijkf, ifluc
291 
292  REAL(RFREAL) :: sgn, ds, sxn, syn, szn, pr, mach, amp, eps, ravgtim
293  REAL(RFREAL) :: fluc(bcdat_inflow_nelm)
294  REAL(RFREAL), POINTER :: cv(:,:), dv(:,:), gv(:,:), vals(:,:), sface(:,:)
295  REAL(RFREAL), POINTER :: tav(:,:)
296 
297 !******************************************************************************
298 
299  CALL registerfunction( region%global,'RFLO_BcondInflow',&
300  'RFLO_ModBoundaryConditions.F90' )
301 
302 ! get dimensions and pointers
303 
304  ifluc = region%global%infloNijk
305  ilev = region%currLevel
306  lbound = patch%lbound
307 
308  CALL rflo_getpatchindices( region,patch,ilev, &
309  ibeg,iend,jbeg,jend,kbeg,kend )
311  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
312  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
313 
314  noff = abs(patch%l1end-patch%l1beg) + 1
315  bctype = patch%bcType
316  distrib = patch%mixt%distrib
317  bcopttype = patch%mixt%switches(bcswi_inflow_type)
318  bcoptfixed = patch%mixt%switches(bcswi_inflow_fixed)
319  bcoptmodel = patch%mixt%switches(bcswi_inflow_model)
320  amp = patch%mixt%amplitude
321  indcp = region%levels(ilev)%mixt%indCp
322  indmol = region%levels(ilev)%mixt%indMol
323  gasmodel = region%mixtInput%gasModel
324 
325  cv => region%levels(ilev)%mixt%cv
326  dv => region%levels(ilev)%mixt%dv
327  gv => region%levels(ilev)%mixt%gv
328  vals => patch%mixt%vals
329 
330 ! initial fluctuations and pointer for turbulent inflow case
331 
332  fluc = 0._rfreal
333  eps = 100._rfreal*epsilon( 1._rfreal )
334 
335 #ifdef STATS
336  IF (bcoptmodel == bcopt_unsteady) THEN
337  ravgtim = 1._rfreal/(region%global%integrTime + eps)
338  tav => region%levels(ilev)%mixt%tav
339 !if(region%global%myProcid==MASTERPROC) write(*,*)'ok1', region%global%integrTime
340  ENDIF
341 #endif
342 
343 ! to take the right face vector and make it point outwards
344 
345  sgn = +1._rfreal
346  inode = 0
347  jnode = 0
348  knode = 0
349  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
350  sgn = -1._rfreal
351  inode = -idir
352  jnode = -jdir
353  knode = -kdir
354  ENDIF
355 
356 ! get the appropriate face vector
357 
358  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
359  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
360  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
361 
362 ! loop over all cells of the patch
363 
364  DO idum=1,region%nDumCells
365  DO k=kbeg,kend
366  DO j=jbeg,jend
367  DO i=ibeg,iend
368  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
369 
370  IF (idum == 1) THEN
371  ijkc = indijk(i,j,k,icoff,ijcoff)
372  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
373  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
374  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
375  sface(zcoord,ijkn)*sface(zcoord,ijkn))
376  sxn = sgn*sface(xcoord,ijkn)/ds
377  syn = sgn*sface(ycoord,ijkn)/ds
378  szn = sgn*sface(zcoord,ijkn)/ds
379 
380  IF (lbound==1 .OR. lbound==2) THEN
381  n1 = j - jbeg
382  n2 = k - kbeg
383  ELSE IF (lbound==3 .OR. lbound==4) THEN
384  n1 = k - kbeg
385  n2 = i - ibeg
386  ELSE IF (lbound==5 .OR. lbound==6) THEN
387  n1 = i - ibeg
388  n2 = j - jbeg
389  ENDIF
390  i2d = distrib * indij(n1,n2,noff)
391 
392  IF (lbound==1) THEN
393  ijkf = indijk(ibeg+ifluc-1,j,k,icoff,ijcoff)
394  ELSEIF (lbound==2) THEN
395  ijkf = indijk(iend-ifluc+1,j,k,icoff,ijcoff)
396  ELSEIF (lbound==3) THEN
397  ijkf = indijk(i,jbeg+ifluc-1,k,icoff,ijcoff)
398  ELSEIF (lbound==4) THEN
399  ijkf = indijk(i,jbeg-ifluc+1,k,icoff,ijcoff)
400  ELSEIF (lbound==5) THEN
401  ijkf = indijk(i,j,kbeg+ifluc-1,icoff,ijcoff)
402  ELSEIF (lbound==6) THEN
403  ijkf = indijk(i,j,kbeg-ifluc+1,icoff,ijcoff)
404  ENDIF
405 
406  IF (gasmodel == gas_model_tcperf) THEN
407 #ifdef STATS
408  IF ((bctype == bc_inflow_veltemp .OR. &
409  bctype == bc_inflow_velpress) .AND. &
410  (region%global%integrTime > eps) .AND. &
411  (region%global%infloNijk < nijk_inflow_init)) THEN
412  CALL rflo_bcondinflowfluc( bcoptmodel, &
413  fluc(bcdat_inflow_u), &
414  fluc(bcdat_inflow_v), &
415  fluc(bcdat_inflow_w), &
416  fluc(bcdat_inflow_t), &
417  fluc(bcdat_inflow_p), &
418  dv(dv_mixt_uvel,ijkf), &
419  dv(dv_mixt_vvel,ijkf), &
420  dv(dv_mixt_wvel,ijkf), &
421  dv(dv_mixt_temp,ijkf), &
422  dv(dv_mixt_pres,ijkf), &
423  tav(2,ijkf)*ravgtim, &
424  tav(3,ijkf)*ravgtim, &
425  tav(4,ijkf)*ravgtim, &
426  tav(5,ijkf)*ravgtim, &
427  tav(6,ijkf)*ravgtim )
428 !if (region%iRegionGlobal==1 .AND. j==17 .AND. k==17) then
429 !write(*,*)'tav',tav(2,ijkF)*ravgtim,tav(3,ijkF)*ravgtim,tav(4,ijkF)*ravgtim,tav(5,ijkF)*ravgtim,tav(6,ijkF)*ravgtim
430 !write(*,*)'dv ',dv(DV_MIXT_UVEL,ijkF),dv(DV_MIXT_VVEL,ijkF),dv(DV_MIXT_WVEL,ijkF),dv(DV_MIXT_TEMP,ijkF),dv(DV_MIXT_PRES,ijkF)
431 !endif
432  ENDIF
433 #endif
434  IF (bctype == bc_inflow_totang) THEN
435 
436  IF (bcopttype == bcopt_subsonic) THEN
437  mach = 0._rfreal
438  ELSE
439  mach = vals(bcdat_inflow_mach,i2d)
440  ENDIF
441  IF (bcoptmodel == bcopt_unsteady) THEN
442  CALL errorstop( region%global,err_val_bcval,&
443  __line__, &
444  'RECYCTURB > 0 not available for BC_INFLOW_TOTANG' )
445  ENDIF
446  CALL bcondinflowperf( bcopttype,bcoptfixed , &
447  vals(bcdat_inflow_ptot ,i2d), &
448  vals(bcdat_inflow_ttot ,i2d), &
449  vals(bcdat_inflow_betah,i2d), &
450  vals(bcdat_inflow_betav,i2d), &
451  mach,sxn,syn,szn , &
452  gv(gv_mixt_cp ,ijkc*indcp ), &
453  gv(gv_mixt_mol ,ijkc*indmol), &
454  cv(cv_mixt_dens,ijkc) , &
455  cv(cv_mixt_xmom,ijkc) , &
456  cv(cv_mixt_ymom,ijkc) , &
457  cv(cv_mixt_zmom,ijkc) , &
458  cv(cv_mixt_dens,ijkd) , &
459  cv(cv_mixt_xmom,ijkd) , &
460  cv(cv_mixt_ymom,ijkd) , &
461  cv(cv_mixt_zmom,ijkd) , &
462  cv(cv_mixt_ener,ijkd),pr )
463 
464  ELSEIF (bctype == bc_inflow_veltemp) THEN
465 
466  CALL rflo_bcondinflowvtperf( bcoptmodel, bcopttype, &
467  vals(bcdat_inflow_u ,i2d)+amp*fluc(bcdat_inflow_u), &
468  vals(bcdat_inflow_v ,i2d)+amp*fluc(bcdat_inflow_v), &
469  vals(bcdat_inflow_w ,i2d)+amp*fluc(bcdat_inflow_w), &
470  vals(bcdat_inflow_t ,i2d)+amp*fluc(bcdat_inflow_t), &
471  vals(bcdat_inflow_p ,i2d)+amp*fluc(bcdat_inflow_p), &
472  sxn,syn,szn , &
473  gv(gv_mixt_cp ,ijkc*indcp),gv(gv_mixt_mol ,ijkc*indmol), &
474  cv(cv_mixt_dens,ijkc) ,cv(cv_mixt_xmom,ijkc) , &
475  cv(cv_mixt_ymom,ijkc) ,cv(cv_mixt_zmom,ijkc) , &
476  cv(cv_mixt_ener,ijkc) ,dv(dv_mixt_pres,ijkc) , &
477  cv(cv_mixt_dens,ijkd) ,cv(cv_mixt_xmom,ijkd) , &
478  cv(cv_mixt_ymom,ijkd) ,cv(cv_mixt_zmom,ijkd) , &
479  cv(cv_mixt_ener,ijkd) )
480 #ifdef STATS
481 ! IF ((region%global%integrTime > eps) .AND. &
482 ! (region%global%infloNijk < NIJK_INFLOW_INIT)) THEN
483 ! cv(CV_MIXT_XMOM,ijkD) = cv(CV_MIXT_XMOM,ijkD) + &
484 ! (vals(BCDAT_INFLOW_U,i2d) - tav(2,ijkD)*rAvgTim)* &
485 ! cv(CV_MIXT_DENS,ijkD)
486 ! ENDIF
487 #endif
488  ELSEIF (bctype == bc_inflow_velpress) THEN
489 
490  CALL rflo_bcondinflowvpperf( bcoptmodel, bcopttype, &
491  vals(bcdat_inflow_u ,i2d)+amp*fluc(bcdat_inflow_u), &
492  vals(bcdat_inflow_v ,i2d)+amp*fluc(bcdat_inflow_v), &
493  vals(bcdat_inflow_w ,i2d)+amp*fluc(bcdat_inflow_w), &
494  vals(bcdat_inflow_t ,i2d)+amp*fluc(bcdat_inflow_t), &
495  vals(bcdat_inflow_p ,i2d)+amp*fluc(bcdat_inflow_p), &
496  sxn,syn,szn , &
497  gv(gv_mixt_cp ,ijkc*indcp),gv(gv_mixt_mol ,ijkc*indmol), &
498  cv(cv_mixt_dens,ijkc) ,cv(cv_mixt_xmom,ijkc) , &
499  cv(cv_mixt_ymom,ijkc) ,cv(cv_mixt_zmom,ijkc) , &
500  cv(cv_mixt_ener,ijkc) ,dv(dv_mixt_pres,ijkc) , &
501  cv(cv_mixt_dens,ijkd) ,cv(cv_mixt_xmom,ijkd) , &
502  cv(cv_mixt_ymom,ijkd) ,cv(cv_mixt_zmom,ijkd) , &
503  cv(cv_mixt_ener,ijkd) )
504 #ifdef STATS
505 ! IF ((region%global%integrTime > eps) .AND. &
506 ! (region%global%infloNijk < NIJK_INFLOW_INIT)) THEN
507 ! cv(CV_MIXT_XMOM,ijkD) = cv(CV_MIXT_XMOM,ijkD) + &
508 ! (vals(BCDAT_INFLOW_U,i2d) - tav(2,ijkD)*rAvgTim)* &
509 ! cv(CV_MIXT_DENS,ijkD)
510 ! ENDIF
511 #endif
512  ENDIF
513  ELSE
514  CALL errorstop( region%global,err_unknown_bc,&
515  __line__ )
516  ENDIF
517 
518  ELSE ! idum > 1
519  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
520  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
521 
522  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
523  cv(cv_mixt_dens,ijkc1)
524  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
525  cv(cv_mixt_xmom,ijkc1)
526  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
527  cv(cv_mixt_ymom,ijkc1)
528  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
529  cv(cv_mixt_zmom,ijkc1)
530  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
531  cv(cv_mixt_ener,ijkc1)
532  ENDIF
533 
534  IF (gasmodel == gas_model_tcperf) THEN
535  CALL mixtureproperties( region,ijkd,ijkd,.false. )
536  ELSE
537  CALL mixtureproperties( region,ijkd,ijkd,.true. )
538  ENDIF
539 
540  ENDDO ! i
541  ENDDO ! j
542  ENDDO ! k
543  ENDDO ! idum
544 
545 ! finalize
546 
547  CALL deregisterfunction( region%global )
548 
549 END SUBROUTINE rflo_bcondinflow
550 
551 !******************************************************************************
552 !
553 ! Purpose: compute mfrate at injection boundary patch by aP^n method.
554 !
555 ! Description: mfrate = solid-density*a*pressure^n
556 !
557 ! Input: region = region dimensions, user input
558 ! patch = current patch.
559 !
560 ! Output: patch%mixt%vals(BCDAT_INJECT_MFRATE,:) = mfrate at current patch
561 !
562 !******************************************************************************
563 
564 SUBROUTINE rflo_bcondinjectionapn( region,patch )
565 
567  IMPLICIT NONE
568 
569 #include "Indexing.h"
570 
571 ! ... parameters
572  TYPE(t_region) :: region
573  TYPE(t_patch) :: patch
574 
575 ! ... loop variables
576  INTEGER :: i, j, k, n1, n2
577 
578 ! ... local variables
579  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, ilev, lbound
580  INTEGER :: idir, jdir, kdir, inode, jnode, knode
581  INTEGER :: icoff, ijcoff, inoff, ijnoff, noff, ijkc, ijkn, i2d
582  INTEGER :: distrib
583 
584  REAL(RFREAL) :: pb, sdens, coeff, power
585  REAL(RFREAL) :: sgn, ds, mrate, sxn, syn, szn
586  REAL(RFREAL), POINTER :: dv(:,:), vals(:,:), sface(:,:)
587 
588 !******************************************************************************
589 
590  CALL registerfunction( region%global,'RFLO_BcondInjectionAPN',&
591  'RFLO_ModBoundaryConditions.F90' )
592 
593 ! get dimensions and pointers -------------------------------------------------
594 
595  ilev = region%currLevel
596  lbound = patch%lbound
597 
598  CALL rflo_getpatchindices( region,patch,ilev, &
599  ibeg,iend,jbeg,jend,kbeg,kend )
601  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
602  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
603 
604  noff = abs(patch%l1end-patch%l1beg) + 1
605  distrib = patch%mixt%distrib
606 
607  dv => region%levels(ilev)%mixt%dv
608  vals => patch%mixt%vals
609 
610 ! to take the right face vector and make it point outwards
611 
612  sgn = +1._rfreal
613  inode = 0
614  jnode = 0
615  knode = 0
616  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
617  sgn = -1._rfreal
618  inode = -idir
619  jnode = -jdir
620  knode = -kdir
621  ENDIF
622 
623 ! get the appropriate face vector
624 
625  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
626  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
627  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
628 
629 ! loop over patch faces -------------------------------------------------------
630 
631  DO k=kbeg,kend
632  DO j=jbeg,jend
633  DO i=ibeg,iend
634  ijkc = indijk(i,j,k,icoff,ijcoff)
635  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
636  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
637  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
638  sface(zcoord,ijkn)*sface(zcoord,ijkn))
639  sxn = sgn*sface(xcoord,ijkn)/ds
640  syn = sgn*sface(ycoord,ijkn)/ds
641  szn = sgn*sface(zcoord,ijkn)/ds
642 
643  IF (lbound==1 .OR. lbound==2) THEN
644  n1 = j - jbeg
645  n2 = k - kbeg
646  ELSE IF (lbound==3 .OR. lbound==4) THEN
647  n1 = k - kbeg
648  n2 = i - ibeg
649  ELSE IF (lbound==5 .OR. lbound==6) THEN
650  n1 = i - ibeg
651  n2 = j - jbeg
652  ENDIF
653  i2d = distrib * indij(n1,n2,noff)
654  sdens = vals(bcdat_inject_sdens ,i2d)
655  coeff = vals(bcdat_inject_acoeff,i2d)
656  power = vals(bcdat_inject_npower,i2d)
657  pb = dv(dv_mixt_pres,ijkc)
658  vals(bcdat_inject_mfrate,i2d) = sdens*coeff*(pb**power)
659 
660  mrate = vals(bcdat_inject_mfrate,i2d)
661  vals(bcdat_inject_rfvfu,i2d) = -mrate*sxn
662  vals(bcdat_inject_rfvfv,i2d) = -mrate*syn
663  vals(bcdat_inject_rfvfw,i2d) = -mrate*szn
664 
665  ENDDO ! i
666  ENDDO ! j
667  ENDDO ! k
668 
669 ! finalize --------------------------------------------------------------------
670 
671  CALL deregisterfunction( region%global )
672 
673 END SUBROUTINE rflo_bcondinjectionapn
674 
675 !******************************************************************************
676 !
677 ! Purpose: update values in dummy cells at injection boundary patch.
678 !
679 ! Description: the boundary condition uses linear extrapolation
680 ! of conservative variables. In addition, a strong boundary
681 ! condition is applied to the convective fluxes.
682 !
683 ! Input: region = region dimensions, user input
684 ! patch = current patch.
685 !
686 ! Output: region%levels%mixt = flow variables in dummy cells.
687 !
688 ! Notes: if the mass flow rate is zero or negative, the corresponding
689 ! cell face is treated like solid wall.
690 !
691 !******************************************************************************
692 
693 SUBROUTINE rflo_bcondinjection( region,patch )
694 
698  IMPLICIT NONE
699 
700 #include "Indexing.h"
701 
702 ! ... parameters
703  TYPE(t_region) :: region
704  TYPE(t_patch) :: patch
705 
706 ! ... loop variables
707  INTEGER :: idum, i, j, k, n1, n2
708 
709 ! ... local variables
710  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir, ilev, lbound
711  INTEGER :: icoff, ijcoff, inoff, ijnoff, noff, ijkc, ijkc1, ijkd, ijkn, i2d
712  INTEGER :: distrib, gasmodel, flowmodel, dumextrapol, bctype
713  INTEGER :: indcp, indmol, inode, jnode, knode
714 
715  REAL(RFREAL) :: sgn, rhoa, rhoua, rhova, rhowa, rhoea, pa
716  REAL(RFREAL) :: mrate, tburn, rhovrel(3), rgas, ds, sxn, syn, szn
717  REAL(RFREAL) :: uinj, vinj, winj, rhodum, rhoedum, maxchange
718  REAL(RFREAL), POINTER :: cv(:,:), dv(:,:), gv(:,:), vals(:,:), sface(:,:)
719 
720 !******************************************************************************
721 
722  CALL registerfunction( region%global,'RFLO_BcondInjection',&
723  'RFLO_ModBoundaryConditions.F90' )
724 
725 ! get dimensions and pointers -------------------------------------------------
726 
727  ilev = region%currLevel
728  lbound = patch%lbound
729 
730  CALL rflo_getpatchindices( region,patch,ilev, &
731  ibeg,iend,jbeg,jend,kbeg,kend )
733  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
734  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
735 
736  noff = abs(patch%l1end-patch%l1beg) + 1
737  bctype = patch%bcType
738  distrib = patch%mixt%distrib
739  dumextrapol = patch%mixt%switches(bcswi_inject_extrap)
740  maxchange = patch%mixt%maxChange
741  gasmodel = region%mixtInput%gasModel
742  flowmodel = region%mixtInput%flowModel
743  indcp = region%levels(ilev)%mixt%indCp
744  indmol = region%levels(ilev)%mixt%indMol
745 
746  cv => region%levels(ilev)%mixt%cv
747  dv => region%levels(ilev)%mixt%dv
748  gv => region%levels(ilev)%mixt%gv
749  vals => patch%mixt%vals
750 
751 ! to take the right face vector and make it point outwards
752 
753  sgn = +1._rfreal
754  inode = 0
755  jnode = 0
756  knode = 0
757  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
758  sgn = -1._rfreal
759  inode = -idir
760  jnode = -jdir
761  knode = -kdir
762  ENDIF
763 
764 ! get the appropriate face vector
765 
766  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
767  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
768  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
769 
770 ! loop over patch and dummy cells ---------------------------------------------
771 
772  DO idum=1,region%nDumCells
773  DO k=kbeg,kend
774  DO j=jbeg,jend
775  DO i=ibeg,iend
776  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
777  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
778  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
779  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
780  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
781  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
782  sface(zcoord,ijkn)*sface(zcoord,ijkn))
783  sxn = sgn*sface(xcoord,ijkn)/ds
784  syn = sgn*sface(ycoord,ijkn)/ds
785  szn = sgn*sface(zcoord,ijkn)/ds
786 
787  IF (lbound==1 .OR. lbound==2) THEN
788  n1 = j - jbeg
789  n2 = k - kbeg
790  ELSE IF (lbound==3 .OR. lbound==4) THEN
791  n1 = k - kbeg
792  n2 = i - ibeg
793  ELSE IF (lbound==5 .OR. lbound==6) THEN
794  n1 = i - ibeg
795  n2 = j - jbeg
796  ENDIF
797  i2d = distrib * indij(n1,n2,noff)
798  mrate = vals(bcdat_inject_mfrate,i2d)
799  tburn = vals(bcdat_inject_temp ,i2d)
800  rhovrel(1) = vals(bcdat_inject_rfvfu ,i2d)
801  rhovrel(2) = vals(bcdat_inject_rfvfv ,i2d)
802  rhovrel(3) = vals(bcdat_inject_rfvfw ,i2d)
803 
804 ! ------- surface burning and 1st dummmy cell
805 ! ------- or heat transfer wall
806 
807  IF ((mrate>0._rfreal .AND. idum==1) .OR. &
808  (bctype == bc_injection_ht)) THEN
809  IF (gasmodel == gas_model_tcperf) THEN
810  CALL bcondinjectionperf( distrib,mrate,tburn,rhovrel,sxn,syn,szn,&
811  gv(gv_mixt_cp ,ijkc*indcp ), &
812  gv(gv_mixt_mol ,ijkc*indmol), &
813  dv(dv_mixt_pres,ijkc ), &
814  rhoa,rhoua,rhova,rhowa,rhoea,pa, &
815  uinj,vinj,winj )
816  ELSE
817  CALL errorstop( region%global,err_unknown_bc,&
818  __line__ )
819  ENDIF
820 
821  rhodum = 2._rfreal*rhoa - cv(cv_mixt_dens,ijkc)
822  rhoedum = 2._rfreal*rhoea - cv(cv_mixt_ener,ijkc)
823  IF (rhodum<=0._rfreal .OR. rhoedum<=0._rfreal .OR. &
824  abs(dv(dv_mixt_temp,ijkc)-tburn)>maxchange*tburn .OR. &
825  dumextrapol==extrapol_const) THEN
826  cv(cv_mixt_dens,ijkd) = rhoa
827  cv(cv_mixt_xmom,ijkd) = rhoua
828  cv(cv_mixt_ymom,ijkd) = rhova
829  cv(cv_mixt_zmom,ijkd) = rhowa
830  cv(cv_mixt_ener,ijkd) = rhoea
831  ELSE
832  cv(cv_mixt_dens,ijkd) = 2._rfreal*rhoa - cv(cv_mixt_dens,ijkc)
833  cv(cv_mixt_xmom,ijkd) = 2._rfreal*rhoua - cv(cv_mixt_xmom,ijkc)
834  cv(cv_mixt_ymom,ijkd) = 2._rfreal*rhova - cv(cv_mixt_ymom,ijkc)
835  cv(cv_mixt_zmom,ijkd) = 2._rfreal*rhowa - cv(cv_mixt_zmom,ijkc)
836  cv(cv_mixt_ener,ijkd) = 2._rfreal*rhoea - cv(cv_mixt_ener,ijkc)
837  ENDIF
838 
839 ! ------- not burning or dummy cell>1, inviscid flow
840 
841  ELSE IF (flowmodel == flow_euler) THEN
842  rhodum = 2._rfreal*cv(cv_mixt_dens,ijkc) - cv(cv_mixt_dens,ijkc1)
843  rhoedum = 2._rfreal*cv(cv_mixt_ener,ijkc) - cv(cv_mixt_ener,ijkc1)
844  IF (rhodum<=0._rfreal .OR. rhoedum<=0._rfreal .OR. &
845  abs(cv(cv_mixt_dens,ijkd)-rhodum)> &
846  maxchange*cv(cv_mixt_dens,ijkd) .OR. &
847  abs(cv(cv_mixt_ener,ijkd)-rhoedum)> &
848  maxchange*cv(cv_mixt_ener,ijkd) .OR. &
849  dumextrapol==extrapol_const) THEN
850  cv(cv_mixt_dens,ijkd) = cv(cv_mixt_dens,ijkc)
851  cv(cv_mixt_xmom,ijkd) = cv(cv_mixt_xmom,ijkc)
852  cv(cv_mixt_ymom,ijkd) = cv(cv_mixt_ymom,ijkc)
853  cv(cv_mixt_zmom,ijkd) = cv(cv_mixt_zmom,ijkc)
854  cv(cv_mixt_ener,ijkd) = cv(cv_mixt_ener,ijkc)
855  ELSE
856  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
857  cv(cv_mixt_dens,ijkc1)
858  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
859  cv(cv_mixt_xmom,ijkc1)
860  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
861  cv(cv_mixt_ymom,ijkc1)
862  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
863  cv(cv_mixt_zmom,ijkc1)
864  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
865  cv(cv_mixt_ener,ijkc1)
866  ENDIF
867 
868 ! ------- not burning or dummy cell>1, viscous flow
869 ! ------- or heat transfer wall
870 
871  ELSE IF (flowmodel == flow_navst) THEN
872  IF ((mrate > 0._rfreal) .OR. &
873  (bctype == bc_injection_ht)) THEN ! burning => extrapolate
874  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
875  cv(cv_mixt_dens,ijkc1)
876  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
877  cv(cv_mixt_xmom,ijkc1)
878  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
879  cv(cv_mixt_ymom,ijkc1)
880  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
881  cv(cv_mixt_zmom,ijkc1)
882  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
883  cv(cv_mixt_ener,ijkc1)
884  ELSE ! not burning => like noslip wall
885  ijkc = indijk(i+(idum-1)*idir,j+(idum-1)*jdir,k+(idum-1)*kdir,icoff,ijcoff)
886  cv(cv_mixt_dens,ijkd) = cv(cv_mixt_dens,ijkc)
887  cv(cv_mixt_xmom,ijkd) = -cv(cv_mixt_xmom,ijkc)
888  cv(cv_mixt_ymom,ijkd) = -cv(cv_mixt_ymom,ijkc)
889  cv(cv_mixt_zmom,ijkd) = -cv(cv_mixt_zmom,ijkc)
890  cv(cv_mixt_ener,ijkd) = cv(cv_mixt_ener,ijkc)
891  ENDIF
892  ENDIF
893 
894  IF (gasmodel == gas_model_tcperf) THEN
895  CALL mixtureproperties( region,ijkd,ijkd,.false. )
896  ELSE
897  CALL mixtureproperties( region,ijkd,ijkd,.true. )
898  ENDIF
899 
900  ENDDO ! i
901  ENDDO ! j
902  ENDDO ! k
903  ENDDO ! idum
904 
905 ! finalize --------------------------------------------------------------------
906 
907  CALL deregisterfunction( region%global )
908 
909 END SUBROUTINE rflo_bcondinjection
910 
911 !******************************************************************************
912 !
913 ! Purpose: convert mdot into density*velocity for injection boundaries.
914 !
915 ! Description: none.
916 !
917 ! Input: .
918 !
919 ! Output: regions = BC data.
920 !
921 ! Notes: none.
922 !
923 !******************************************************************************
924 
925 SUBROUTINE rflo_bcondinjectioninit( region )
926 
929  IMPLICIT NONE
930 
931 #include "Indexing.h"
932 
933 ! ... parameters
934  TYPE(t_region) :: region
935 
936 ! ... loop variables
937  INTEGER :: ilev, ipatch, i, j, k, n1, n2
938 
939 ! ... local variables
940  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir, lbound
941  INTEGER :: icoff, ijcoff, inoff, ijnoff, noff, ijkc, ijkn, i2d
942  INTEGER :: inode, jnode, knode
943 
944  REAL(RFREAL) :: sgn, ds, sxn, syn, szn, mrate
945  REAL(RFREAL) :: sdens, coeff, power, pb
946  REAL(RFREAL), POINTER :: dv(:,:), vals(:,:), sface(:,:)
947 
948  TYPE(t_patch), POINTER :: patch, patchprev
949  TYPE(t_global), POINTER :: global
950 
951 !******************************************************************************
952 
953  global => region%global
954 
955  CALL registerfunction( global,'RFLO_BcondInjectionInit',&
956  'RFLO_ModBoundaryConditions.F90' )
957 
958 ! copy values/distribution to variables ---------------------------------------
959 
960  DO ilev=1,region%nGridLevels
961  dv => region%levels(ilev)%mixt%dv
962 
963  DO ipatch=1,region%nPatches
964 
965  patch => region%levels(ilev)%patches(ipatch)
966  IF (ilev > 1) patchprev => region%levels(ilev-1)%patches(ipatch)
967 
968  IF (patch%bcType>=bc_injection .AND. &
969  patch%bcType<=bc_injection+bc_range) THEN
970 
971  lbound = patch%lbound
972  noff = abs(patch%l1end-patch%l1beg) + 1
973  vals => patch%mixt%vals
974 
975 ! ----- BC values as distribution
976 
977  IF (patch%mixt%distrib==bcdat_distrib) THEN
978  CALL rflo_getpatchindices( region,patch,ilev, &
979  ibeg,iend,jbeg,jend,kbeg,kend )
981  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
982  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
983 
984 ! ------- to take the right face vector and make it point outwards
985 
986  sgn = +1._rfreal
987  inode = 0
988  jnode = 0
989  knode = 0
990  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
991  sgn = -1._rfreal
992  inode = -idir
993  jnode = -jdir
994  knode = -kdir
995  ENDIF
996 
997 ! ------- get the appropriate face vector
998 
999  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
1000  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
1001  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
1002 
1003 ! ------- mRate for injection APN
1004 
1005  IF (patch%bcType==bc_injection_apn) THEN
1006  DO k=kbeg,kend
1007  DO j=jbeg,jend
1008  DO i=ibeg,iend
1009  ijkc = indijk(i,j,k,icoff,ijcoff)
1010 
1011  IF (lbound==1 .OR. lbound==2) THEN
1012  n1 = j - jbeg
1013  n2 = k - kbeg
1014  ELSE IF (lbound==3 .OR. lbound==4) THEN
1015  n1 = k - kbeg
1016  n2 = i - ibeg
1017  ELSE IF (lbound==5 .OR. lbound==6) THEN
1018  n1 = i - ibeg
1019  n2 = j - jbeg
1020  ENDIF
1021  i2d = indij(n1,n2,noff)
1022  sdens = vals(bcdat_inject_sdens ,i2d)
1023  coeff = vals(bcdat_inject_acoeff,i2d)
1024  power = vals(bcdat_inject_npower,i2d)
1025  pb = dv(dv_mixt_pres,ijkc)
1026  vals(bcdat_inject_mfrate,i2d) = sdens*coeff*(pb**power)
1027 
1028  ENDDO ! i
1029  ENDDO ! j
1030  ENDDO ! k
1031  ENDIF ! bcType
1032 
1033 ! ------- loop over the cells
1034 
1035  DO k=kbeg,kend
1036  DO j=jbeg,jend
1037  DO i=ibeg,iend
1038  ijkc = indijk(i,j,k,icoff,ijcoff)
1039  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
1040  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
1041  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
1042  sface(zcoord,ijkn)*sface(zcoord,ijkn))
1043  sxn = sgn*sface(xcoord,ijkn)/ds
1044  syn = sgn*sface(ycoord,ijkn)/ds
1045  szn = sgn*sface(zcoord,ijkn)/ds
1046  IF (lbound==1 .OR. lbound==2) THEN
1047  n1 = j - jbeg
1048  n2 = k - kbeg
1049  ELSE IF (lbound==3 .OR. lbound==4) THEN
1050  n1 = k - kbeg
1051  n2 = i - ibeg
1052  ELSE IF (lbound==5 .OR. lbound==6) THEN
1053  n1 = i - ibeg
1054  n2 = j - jbeg
1055  ENDIF
1056  i2d = indij(n1,n2,noff)
1057  mrate = vals(bcdat_inject_mfrate,i2d)
1058  vals(bcdat_inject_rfvfu,i2d) = -mrate*sxn
1059  vals(bcdat_inject_rfvfv,i2d) = -mrate*syn
1060  vals(bcdat_inject_rfvfw,i2d) = -mrate*szn
1061  ENDDO ! i
1062  ENDDO ! j
1063  ENDDO ! k
1064 
1065 ! ----- BC values constant
1066 
1067  ELSE
1068  vals(bcdat_inject_rfvfu,:) = 0._rfreal
1069  vals(bcdat_inject_rfvfv,:) = 0._rfreal
1070  vals(bcdat_inject_rfvfw,:) = 0._rfreal
1071  ENDIF ! distrib
1072 
1073  IF (ilev > 1) THEN
1074  CALL rflo_copyboundarydata( global,patchprev,patch )
1075  ENDIF
1076 
1077  ENDIF ! bcType
1078 
1079  ENDDO ! iPatch
1080  ENDDO ! iLev
1081 
1082 ! finalize --------------------------------------------------------------------
1083 
1084  CALL deregisterfunction( global )
1085 
1086 END SUBROUTINE rflo_bcondinjectioninit
1087 
1088 !******************************************************************************
1089 !
1090 ! Purpose: update values in dummy cells at noslip (viscous) walls.
1091 !
1092 ! Description: none.
1093 !
1094 ! Input: region = region dimensions, user input
1095 ! patch = current patch.
1096 !
1097 ! Output: region%levels%mixt = flow variables in dummy cells.
1098 !
1099 ! Notes: none.
1100 !
1101 !******************************************************************************
1102 
1103 SUBROUTINE rflo_bcondnoslipwall( region,patch )
1104 
1108  IMPLICIT NONE
1109 
1110 #include "Indexing.h"
1111 
1112 ! ... parameters
1113  TYPE(t_region) :: region
1114  TYPE(t_patch) :: patch
1115 
1116 ! ... loop variables
1117  INTEGER :: idum, i, j, k, n1, n2
1118 
1119 ! ... local variables
1120  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
1121  INTEGER :: ilev, lbound, icoff, ijcoff, noff, distrib, bcopt
1122  INTEGER :: indcp, indmol, gasmodel, ijkc, ijkc1, ijkci, ijkcb, ijkd, i2d
1123  INTEGER :: inode, jnode, knode, inoff, ijnoff, ijkn(4)
1124 
1125  LOGICAL :: movegrid
1126 
1127  REAL(RFREAL) :: twall, rgas, gamma, pb, u, v, w, temp, rho, rhoe
1128  REAL(RFREAL) :: dt, dxn(4), dyn(4), dzn(4), usurf, vsurf, wsurf
1129  REAL(RFREAL), POINTER :: cv(:,:), dv(:,:), gv(:,:), vals(:,:)
1130  REAL(RFREAL), POINTER :: xyz(:,:), xyzold(:,:)
1131 
1132 !******************************************************************************
1133 
1134  CALL registerfunction( region%global,'RFLO_BcondNoslipWall',&
1135  'RFLO_ModBoundaryConditions.F90' )
1136 
1137 ! get dimensions and pointers -------------------------------------------------
1138 
1139  ilev = region%currLevel
1140  lbound = patch%lbound
1141 
1142  CALL rflo_getpatchindices( region,patch,ilev, &
1143  ibeg,iend,jbeg,jend,kbeg,kend )
1145  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
1146  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
1147 
1148  noff = abs(patch%l1end-patch%l1beg) + 1
1149  distrib = patch%mixt%distrib
1150  bcopt = patch%mixt%switches(bcswi_noslip_adiabat)
1151  indcp = region%levels(ilev)%mixt%indCp
1152  indmol = region%levels(ilev)%mixt%indMol
1153  gasmodel = region%mixtInput%gasModel
1154  movegrid = region%mixtInput%moveGrid
1155  dt = region%global%dtMin
1156 
1157  cv => region%levels(ilev)%mixt%cv
1158  dv => region%levels(ilev)%mixt%dv
1159  gv => region%levels(ilev)%mixt%gv
1160  xyz => region%levels(ilev)%grid%xyz
1161  IF (movegrid) xyzold => region%levels(ilev)%gridOld%xyz
1162  vals => patch%mixt%vals
1163 
1164 ! to take the correct surface coordinates
1165 
1166  inode = 0
1167  jnode = 0
1168  knode = 0
1169  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
1170  inode = -idir
1171  jnode = -jdir
1172  knode = -kdir
1173  ENDIF
1174 
1175 ! loop over all cells of a patch ----------------------------------------------
1176 
1177  DO idum=1,region%nDumCells
1178  DO k=kbeg,kend
1179  DO j=jbeg,jend
1180  DO i=ibeg,iend
1181  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
1182  ijkci = indijk(i+(idum-1)*idir,j+(idum-1)*jdir,k+(idum-1)*kdir,icoff,ijcoff)
1183 
1184 ! ------- wall velocity in case of moving surface
1185 
1186  IF (movegrid) THEN
1187  ijkn(1) = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
1188  IF (lbound==1 .OR. lbound==2) THEN
1189  ijkn(2) = indijk(i+inode,j+jnode+1,k+knode ,inoff,ijnoff)
1190  ijkn(3) = indijk(i+inode,j+jnode+1,k+knode+1,inoff,ijnoff)
1191  ijkn(4) = indijk(i+inode,j+jnode ,k+knode+1,inoff,ijnoff)
1192  ELSE IF (lbound==3 .OR. lbound==4) THEN
1193  ijkn(2) = indijk(i+inode ,j+jnode,k+knode+1,inoff,ijnoff)
1194  ijkn(3) = indijk(i+inode+1,j+jnode,k+knode+1,inoff,ijnoff)
1195  ijkn(4) = indijk(i+inode+1,j+jnode,k+knode ,inoff,ijnoff)
1196  ELSE IF (lbound==5 .OR. lbound==6) THEN
1197  ijkn(2) = indijk(i+inode+1,j+jnode ,k+knode,inoff,ijnoff)
1198  ijkn(3) = indijk(i+inode+1,j+jnode+1,k+knode,inoff,ijnoff)
1199  ijkn(4) = indijk(i+inode ,j+jnode+1,k+knode,inoff,ijnoff)
1200  ENDIF
1201  dxn(1) = xyz(xcoord,ijkn(1)) - xyzold(xcoord,ijkn(1))
1202  dyn(1) = xyz(ycoord,ijkn(1)) - xyzold(ycoord,ijkn(1))
1203  dzn(1) = xyz(zcoord,ijkn(1)) - xyzold(zcoord,ijkn(1))
1204  dxn(2) = xyz(xcoord,ijkn(2)) - xyzold(xcoord,ijkn(2))
1205  dyn(2) = xyz(ycoord,ijkn(2)) - xyzold(ycoord,ijkn(2))
1206  dzn(2) = xyz(zcoord,ijkn(2)) - xyzold(zcoord,ijkn(2))
1207  dxn(3) = xyz(xcoord,ijkn(3)) - xyzold(xcoord,ijkn(3))
1208  dyn(3) = xyz(ycoord,ijkn(3)) - xyzold(ycoord,ijkn(3))
1209  dzn(3) = xyz(zcoord,ijkn(3)) - xyzold(zcoord,ijkn(3))
1210  dxn(4) = xyz(xcoord,ijkn(4)) - xyzold(xcoord,ijkn(4))
1211  dyn(4) = xyz(ycoord,ijkn(4)) - xyzold(ycoord,ijkn(4))
1212  dzn(4) = xyz(zcoord,ijkn(4)) - xyzold(zcoord,ijkn(4))
1213  usurf = 0.25_rfreal*(dxn(1)+dxn(2)+dxn(3)+dxn(4))/dt
1214  vsurf = 0.25_rfreal*(dyn(1)+dyn(2)+dyn(3)+dyn(4))/dt
1215  wsurf = 0.25_rfreal*(dzn(1)+dzn(2)+dzn(3)+dzn(4))/dt
1216  ELSE
1217  usurf = 0._rfreal
1218  vsurf = 0._rfreal
1219  wsurf = 0._rfreal
1220  ENDIF
1221 
1222 ! ------- adiabatic wall
1223 
1224  IF (bcopt == bcopt_adiabat) THEN
1225  cv(cv_mixt_dens,ijkd) = cv(cv_mixt_dens,ijkci)
1226  cv(cv_mixt_xmom,ijkd) = 2._rfreal*usurf - cv(cv_mixt_xmom,ijkci)
1227  cv(cv_mixt_ymom,ijkd) = 2._rfreal*vsurf - cv(cv_mixt_ymom,ijkci)
1228  cv(cv_mixt_zmom,ijkd) = 2._rfreal*wsurf - cv(cv_mixt_zmom,ijkci)
1229  cv(cv_mixt_ener,ijkd) = cv(cv_mixt_ener,ijkci)
1230 
1231 ! ------- prescribed wall temperature
1232 
1233  ELSE
1234  ijkcb = indijk(i,j,k,icoff,ijcoff)
1235  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
1236  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
1237 
1238  IF (lbound==1 .OR. lbound==2) THEN
1239  n1 = j - jbeg
1240  n2 = k - kbeg
1241  ELSE IF (lbound==3 .OR. lbound==4) THEN
1242  n1 = k - kbeg
1243  n2 = i - ibeg
1244  ELSE IF (lbound==5 .OR. lbound==6) THEN
1245  n1 = i - ibeg
1246  n2 = j - jbeg
1247  ENDIF
1248  i2d = distrib * indij(n1,n2,noff)
1249 
1250  twall = vals(bcdat_noslip_twall,i2d)
1251  pb = dv(dv_mixt_pres,ijkcb)
1252  u = cv(cv_mixt_xmom,ijkci)/cv(cv_mixt_dens,ijkci)
1253  v = cv(cv_mixt_ymom,ijkci)/cv(cv_mixt_dens,ijkci)
1254  w = cv(cv_mixt_zmom,ijkci)/cv(cv_mixt_dens,ijkci)
1255  IF (idum == 1) THEN
1256  temp = 2._rfreal*twall - dv(dv_mixt_temp,ijkc)
1257  ELSE
1258  temp = 2._rfreal*dv(dv_mixt_temp,ijkc) - dv(dv_mixt_temp,ijkc1)
1259  ENDIF
1260 
1261  IF (gasmodel == gas_model_tcperf) THEN
1262  rgas = mixtperf_r_m( gv(gv_mixt_mol,ijkcb*indmol) )
1263  gamma = mixtperf_g_cpr( gv(gv_mixt_cp,ijkcb*indcp),rgas )
1264  rho = mixtperf_d_prt( pb,rgas,temp )
1265  rhoe = rho * mixtperf_eo_dgpuvw( rho,gamma,pb,u,v,w )
1266  ELSE
1267  CALL errorstop( region%global,err_unknown_bc,&
1268  __line__ )
1269  ENDIF
1270  cv(cv_mixt_dens,ijkd) = rho
1271  cv(cv_mixt_xmom,ijkd) = 2._rfreal*usurf - u*rho
1272  cv(cv_mixt_ymom,ijkd) = 2._rfreal*vsurf - v*rho
1273  cv(cv_mixt_zmom,ijkd) = 2._rfreal*wsurf - w*rho
1274  cv(cv_mixt_ener,ijkd) = rhoe
1275  ENDIF ! bcOpt
1276 
1277  IF (gasmodel == gas_model_tcperf) THEN
1278  CALL mixtureproperties( region,ijkd,ijkd,.false. )
1279  ELSE
1280  CALL mixtureproperties( region,ijkd,ijkd,.true. )
1281  ENDIF
1282 
1283  ENDDO ! i
1284  ENDDO ! j
1285  ENDDO ! k
1286  ENDDO ! idum
1287 
1288 ! finalize --------------------------------------------------------------------
1289 
1290  CALL deregisterfunction( region%global )
1291 
1292 END SUBROUTINE rflo_bcondnoslipwall
1293 
1294 !******************************************************************************
1295 !
1296 ! Purpose: update values in dummy cells at outflow boundary patch.
1297 !
1298 ! Description: none.
1299 !
1300 ! Input: region = region dimensions, user input
1301 ! patch = current patch.
1302 !
1303 ! Output: region%levels%mixt = flow variables in dummy cells.
1304 !
1305 ! Notes: none.
1306 !
1307 !******************************************************************************
1308 
1309 SUBROUTINE rflo_bcondoutflow( region,patch )
1310 
1314  IMPLICIT NONE
1315 
1316 #include "Indexing.h"
1317 
1318 ! ... parameters
1319  TYPE(t_region) :: region
1320  TYPE(t_patch) :: patch
1321 
1322 ! ... loop variables
1323  INTEGER :: idum, i, j, k, n1, n2
1324 
1325 ! ... local variables
1326  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
1327  INTEGER :: ilev, lbound, icoff, ijcoff, inoff, ijnoff, noff, distrib, bcopt
1328  INTEGER :: indcp, indmol, gasmodel, ijkc, ijkc1, ijkd, ijkn, i2d
1329  INTEGER :: inode, jnode, knode
1330  INTEGER :: bcmodel
1331 
1332  REAL(RFREAL) :: sgn, ds, sxn, syn, szn, pout
1333  REAL(RFREAL) :: ud, vd, wd, tempd, rold, rrold, uold, vold, wold, &
1334  rdold, rrdold, udold, vdold, wdold, edold, rlen, dtmin, kappa
1335  REAL(RFREAL), POINTER :: cv(:,:), dv(:,:), gv(:,:), vals(:,:), sface(:,:)
1336  REAL(RFREAL), POINTER :: cvold(:,:)
1337 
1338 !******************************************************************************
1339 
1340  CALL registerfunction( region%global,'RFLO_BcondOutflow',&
1341  'RFLO_ModBoundaryConditions.F90' )
1342 
1343 ! get dimensions and pointers
1344 
1345  ilev = region%currLevel
1346  lbound = patch%lbound
1347 
1348  CALL rflo_getpatchindices( region,patch,ilev, &
1349  ibeg,iend,jbeg,jend,kbeg,kend )
1351  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
1352  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
1353 
1354  noff = abs(patch%l1end-patch%l1beg) + 1
1355  distrib = patch%mixt%distrib
1356  bcopt = patch%mixt%switches(bcswi_outflow_type)
1357  bcmodel = patch%mixt%switches(bcswi_outflow_model)
1358  indcp = region%levels(ilev)%mixt%indCp
1359  indmol = region%levels(ilev)%mixt%indMol
1360  gasmodel = region%mixtInput%gasModel
1361 
1362  cv => region%levels(ilev)%mixt%cv
1363  cvold => region%levels(ilev)%mixt%cvOld
1364  dv => region%levels(ilev)%mixt%dv
1365  gv => region%levels(ilev)%mixt%gv
1366  vals => patch%mixt%vals
1367 
1368 ! to take the right face vector and make it point outwards
1369 
1370  sgn = +1._rfreal
1371  inode = 0
1372  jnode = 0
1373  knode = 0
1374  IF (lbound==2 .OR. lbound==4 .OR. lbound==6) THEN
1375  sgn = -1._rfreal
1376  inode = -idir
1377  jnode = -jdir
1378  knode = -kdir
1379  ENDIF
1380 
1381 ! get the appropriate face vector
1382 
1383  IF (lbound==1 .OR. lbound==2) sface => region%levels(ilev)%grid%si
1384  IF (lbound==3 .OR. lbound==4) sface => region%levels(ilev)%grid%sj
1385  IF (lbound==5 .OR. lbound==6) sface => region%levels(ilev)%grid%sk
1386 
1387 ! loop over all cells of a patch
1388 
1389  DO idum=1,region%nDumCells
1390  DO k=kbeg,kend
1391  DO j=jbeg,jend
1392  DO i=ibeg,iend
1393  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
1394 
1395  IF (idum == 1) THEN
1396  ijkc = indijk(i,j,k,icoff,ijcoff)
1397  ijkn = indijk(i+inode,j+jnode,k+knode,inoff,ijnoff)
1398  ds = sqrt(sface(xcoord,ijkn)*sface(xcoord,ijkn)+ &
1399  sface(ycoord,ijkn)*sface(ycoord,ijkn)+ &
1400  sface(zcoord,ijkn)*sface(zcoord,ijkn))
1401  sxn = sgn*sface(xcoord,ijkn)/ds
1402  syn = sgn*sface(ycoord,ijkn)/ds
1403  szn = sgn*sface(zcoord,ijkn)/ds
1404 
1405  IF (lbound==1 .OR. lbound==2) THEN
1406  n1 = j - jbeg
1407  n2 = k - kbeg
1408  ELSE IF (lbound==3 .OR. lbound==4) THEN
1409  n1 = k - kbeg
1410  n2 = i - ibeg
1411  ELSE IF (lbound==5 .OR. lbound==6) THEN
1412  n1 = i - ibeg
1413  n2 = j - jbeg
1414  ENDIF
1415  i2d = distrib * indij(n1,n2,noff)
1416 
1417  IF (gasmodel == gas_model_tcperf) THEN
1418  IF (bcopt == bcopt_supersonic) THEN ! cannot specify pressure
1419  pout = 0._rfreal
1420  ELSE ! take given pressure
1421  pout = vals(bcdat_outflow_press,i2d)
1422  ENDIF
1423 
1424  IF (bcmodel == bcopt_default) THEN
1425  CALL bcondoutflowperf( bcopt,pout,sxn,syn,szn, &
1426  gv(gv_mixt_cp ,ijkc*indcp),gv(gv_mixt_mol ,ijkc*indmol), &
1427  cv(cv_mixt_dens,ijkc) ,cv(cv_mixt_xmom,ijkc) , &
1428  cv(cv_mixt_ymom,ijkc) ,cv(cv_mixt_zmom,ijkc) , &
1429  cv(cv_mixt_ener,ijkc) ,dv(dv_mixt_pres,ijkc) , &
1430  cv(cv_mixt_dens,ijkd) ,cv(cv_mixt_xmom,ijkd) , &
1431  cv(cv_mixt_ymom,ijkd) ,cv(cv_mixt_zmom,ijkd) , &
1432  cv(cv_mixt_ener,ijkd) )
1433  ELSE
1434  ijkc1 = &
1435  indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
1436 
1437  rold = cvold(cv_mixt_dens,ijkc )
1438  rrold = 1._rfreal/rold
1439  uold = cvold(cv_mixt_xmom,ijkc )*rrold
1440  vold = cvold(cv_mixt_ymom,ijkc )*rrold
1441  wold = cvold(cv_mixt_zmom,ijkc )*rrold
1442 
1443  rdold = cvold(cv_mixt_dens,ijkd )
1444  rrdold= 1._rfreal/rdold
1445  udold = cvold(cv_mixt_xmom,ijkd )*rrdold
1446  vdold = cvold(cv_mixt_ymom,ijkd )*rrdold
1447  wdold = cvold(cv_mixt_zmom,ijkd )*rrdold
1448  edold = cvold(cv_mixt_ener,ijkd )*rrdold
1449  rlen = region%global%refLength
1450  dtmin = region%global%dtMin
1451  IF (bcopt == bcopt_supersonic) THEN ! cannot specify nrcoef
1452  kappa = 0._rfreal
1453  ELSE
1454  kappa = vals(bcdat_outflow_nrcoef,i2d)
1455  ENDIF
1456 
1457  ud = 2._rfreal*dv(dv_mixt_uvel,ijkc ) - dv(dv_mixt_uvel,ijkc1)
1458  vd = 2._rfreal*dv(dv_mixt_vvel,ijkc ) - dv(dv_mixt_vvel,ijkc1)
1459  wd = 2._rfreal*dv(dv_mixt_wvel,ijkc ) - dv(dv_mixt_wvel,ijkc1)
1460  tempd = 2._rfreal*dv(dv_mixt_temp,ijkc ) - dv(dv_mixt_temp,ijkc1)
1461 
1462  CALL rflo_bcondoutflowperf( bcopt,pout,sxn,syn,szn, &
1463  gv(gv_mixt_cp ,ijkc*indcp),gv(gv_mixt_mol ,ijkc*indmol), &
1464  cv(cv_mixt_dens,ijkc) ,cv(cv_mixt_xmom,ijkc) , &
1465  cv(cv_mixt_ymom,ijkc) ,cv(cv_mixt_zmom,ijkc) , &
1466  cv(cv_mixt_ener,ijkc) ,dv(dv_mixt_pres,ijkc) , &
1467  ud, vd, wd, tempd, uold, vold, wold , &
1468  udold, vdold, wdold, rdold, edold, rlen, dtmin, kappa , &
1469  cv(cv_mixt_dens,ijkd) ,cv(cv_mixt_xmom,ijkd) , &
1470  cv(cv_mixt_ymom,ijkd) ,cv(cv_mixt_zmom,ijkd) , &
1471  cv(cv_mixt_ener,ijkd) )
1472  ENDIF ! bcModel
1473  ELSE
1474  CALL errorstop( region%global,err_unknown_bc,&
1475  __line__ )
1476  ENDIF ! gasModel
1477 
1478  ELSE ! idum > 1
1479  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
1480  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
1481 
1482  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
1483  cv(cv_mixt_dens,ijkc1)
1484  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
1485  cv(cv_mixt_xmom,ijkc1)
1486  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
1487  cv(cv_mixt_ymom,ijkc1)
1488  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
1489  cv(cv_mixt_zmom,ijkc1)
1490  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
1491  cv(cv_mixt_ener,ijkc1)
1492  ENDIF
1493 
1494  IF (gasmodel == gas_model_tcperf) THEN
1495  CALL mixtureproperties( region,ijkd,ijkd,.false. )
1496  ELSE
1497  CALL mixtureproperties( region,ijkd,ijkd,.true. )
1498  ENDIF
1499 
1500  ENDDO ! i
1501  ENDDO ! j
1502  ENDDO ! k
1503  ENDDO ! idum
1504 
1505 ! finalize
1506 
1507  CALL deregisterfunction( region%global )
1508 
1509 END SUBROUTINE rflo_bcondoutflow
1510 
1511 !******************************************************************************
1512 !
1513 ! Purpose: exchange values between interior and dummy cells of the
1514 ! corresponding patches of a rotationally periodic boundary.
1515 !
1516 ! Description: none.
1517 !
1518 ! Input: region = current region
1519 ! regionSrc = source region
1520 ! patch = current patch of region
1521 ! patchSrc = source patch of regionSrc.
1522 !
1523 ! Output: region%levels%mixt = flow variables in dummy cells.
1524 !
1525 ! Notes: none.
1526 !
1527 !******************************************************************************
1528 
1529 SUBROUTINE rflo_bcondrotatperiod( region,regionSrc,patch,patchSrc )
1530 
1531  USE modinterfaces, ONLY : mixtureproperties
1532  IMPLICIT NONE
1533 
1534 ! ... parameters
1535  TYPE(t_region) :: region, regionsrc
1536  TYPE(t_patch) :: patch, patchsrc
1537 
1538 ! ... loop variables
1539 
1540 
1541 ! ... local variables
1542 
1543 !******************************************************************************
1544 
1545  CALL registerfunction( region%global,'RFLO_BcondRotatPeriod',&
1546  'RFLO_ModBoundaryConditions.F90' )
1547 
1548 
1549 
1550  CALL deregisterfunction( region%global )
1551 
1552 END SUBROUTINE rflo_bcondrotatperiod
1553 
1554 !******************************************************************************
1555 !
1556 ! Purpose: update values in dummy cells at slip walls.
1557 !
1558 ! Description: the boundary condition uses extrapolation of conservative
1559 ! variables (constant or linear). In addition, a strong
1560 ! boundary condition is applied to the convective fluxes.
1561 !
1562 ! Input: region = region dimensions, user input
1563 ! patch = current patch.
1564 !
1565 ! Output: region%levels%mixt = flow variables in dummy cells.
1566 !
1567 ! Notes: none.
1568 !
1569 !******************************************************************************
1570 
1571 SUBROUTINE rflo_bcondslipwall( region,patch )
1572 
1575  IMPLICIT NONE
1576 
1577 #include "Indexing.h"
1578 
1579 ! ... parameters
1580  TYPE(t_region) :: region
1581  TYPE(t_patch) :: patch
1582 
1583 ! ... loop variables
1584  INTEGER :: idum, i, j, k
1585 
1586 ! ... local variables
1587  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
1588  INTEGER :: ilev, icoff, ijcoff, ijkc, ijkc1, ijkd
1589  INTEGER :: dumextrapol, gasmodel
1590 
1591  REAL(RFREAL) :: rhodum, rhoedum, maxchange
1592  REAL(RFREAL) :: uvel, vvel, wvel
1593  REAL(RFREAL), POINTER :: cv(:,:)
1594 
1595 !******************************************************************************
1596 
1597  CALL registerfunction( region%global,'RFLO_BcondSlipWall',&
1598  'RFLO_ModBoundaryConditions.F90' )
1599 
1600 ! get dimensions and pointers
1601 
1602  ilev = region%currLevel
1603 
1604  CALL rflo_getpatchindices( region,patch,ilev, &
1605  ibeg,iend,jbeg,jend,kbeg,kend )
1607  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
1608 
1609  dumextrapol = patch%mixt%switches(bcswi_slipw_extrap)
1610  maxchange = patch%mixt%maxChange
1611  gasmodel = region%mixtInput%gasModel
1612 
1613  cv => region%levels(ilev)%mixt%cv
1614 
1615 ! loop over all cells of a patch
1616 
1617  DO idum=1,region%nDumCells
1618  DO k=kbeg,kend
1619  DO j=jbeg,jend
1620  DO i=ibeg,iend
1621  ijkc = indijk(i-(idum-1)*idir,j-(idum-1)*jdir,k-(idum-1)*kdir,icoff,ijcoff)
1622  ijkc1 = indijk(i-(idum-2)*idir,j-(idum-2)*jdir,k-(idum-2)*kdir,icoff,ijcoff)
1623  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
1624 
1625  rhodum = 2._rfreal*cv(cv_mixt_dens,ijkc) - cv(cv_mixt_dens,ijkc1)
1626  rhoedum = 2._rfreal*cv(cv_mixt_ener,ijkc) - cv(cv_mixt_ener,ijkc1)
1627  IF (rhodum<=0._rfreal .OR. rhoedum<=0._rfreal .OR. &
1628  abs(cv(cv_mixt_dens,ijkd)-rhodum)> &
1629  maxchange*cv(cv_mixt_dens,ijkd) .OR. &
1630  abs(cv(cv_mixt_ener,ijkd)-rhoedum)> &
1631  maxchange*cv(cv_mixt_ener,ijkd) .OR. &
1632  dumextrapol==extrapol_const) THEN
1633  cv(cv_mixt_dens,ijkd) = cv(cv_mixt_dens,ijkc)
1634  cv(cv_mixt_xmom,ijkd) = cv(cv_mixt_xmom,ijkc)
1635  cv(cv_mixt_ymom,ijkd) = cv(cv_mixt_ymom,ijkc)
1636  cv(cv_mixt_zmom,ijkd) = cv(cv_mixt_zmom,ijkc)
1637  cv(cv_mixt_ener,ijkd) = cv(cv_mixt_ener,ijkc)
1638  ELSE
1639  cv(cv_mixt_dens,ijkd) = 2._rfreal*cv(cv_mixt_dens,ijkc ) - &
1640  cv(cv_mixt_dens,ijkc1)
1641  cv(cv_mixt_xmom,ijkd) = 2._rfreal*cv(cv_mixt_xmom,ijkc ) - &
1642  cv(cv_mixt_xmom,ijkc1)
1643  cv(cv_mixt_ymom,ijkd) = 2._rfreal*cv(cv_mixt_ymom,ijkc ) - &
1644  cv(cv_mixt_ymom,ijkc1)
1645  cv(cv_mixt_zmom,ijkd) = 2._rfreal*cv(cv_mixt_zmom,ijkc ) - &
1646  cv(cv_mixt_zmom,ijkc1)
1647  cv(cv_mixt_ener,ijkd) = 2._rfreal*cv(cv_mixt_ener,ijkc ) - &
1648  cv(cv_mixt_ener,ijkc1)
1649  IF (patch%mixt%setMotion .AND. idum==1) THEN
1650  uvel = 2._rfreal*patch%mixt%bndVel(xcoord) - &
1651  cv(cv_mixt_xmom,ijkc)/cv(cv_mixt_dens,ijkc)
1652  vvel = 2._rfreal*patch%mixt%bndVel(ycoord) - &
1653  cv(cv_mixt_ymom,ijkc)/cv(cv_mixt_dens,ijkc)
1654  wvel = 2._rfreal*patch%mixt%bndVel(zcoord) - &
1655  cv(cv_mixt_zmom,ijkc)/cv(cv_mixt_dens,ijkc)
1656  cv(cv_mixt_xmom,ijkd) = uvel*cv(cv_mixt_dens,ijkd )
1657  cv(cv_mixt_ymom,ijkd) = vvel*cv(cv_mixt_dens,ijkd )
1658  cv(cv_mixt_zmom,ijkd) = wvel*cv(cv_mixt_dens,ijkd )
1659  ENDIF
1660  ENDIF
1661 
1662  IF (gasmodel == gas_model_tcperf) THEN
1663  CALL mixtureproperties( region,ijkd,ijkd,.false. )
1664  ELSE
1665  CALL mixtureproperties( region,ijkd,ijkd,.true. )
1666  ENDIF
1667  ENDDO ! i
1668  ENDDO ! j
1669  ENDDO ! k
1670  ENDDO ! idum
1671 
1672 ! finalize
1673 
1674  CALL deregisterfunction( region%global )
1675 
1676 END SUBROUTINE rflo_bcondslipwall
1677 
1678 !******************************************************************************
1679 !
1680 ! Purpose: update values in dummy cells at symmetry boundary patch.
1681 !
1682 ! Description: (1) density and energy in dummy cells are set equal to values
1683 ! in the corresponding interior nodes
1684 !
1685 ! (2) velocity components in dummy cells are mirrored with respect
1686 ! to the boundary; it is assumed that the boundary coincides
1687 ! either with a x=const., y=const., or z=const. plane.
1688 !
1689 ! Input: region = region dimensions, user input
1690 ! patch = current patch.
1691 !
1692 ! Output: region%levels%mixt = flow variables in dummy cells.
1693 !
1694 ! Notes: none.
1695 !
1696 !******************************************************************************
1697 
1698 SUBROUTINE rflo_bcondsymmetry( region,patch )
1699 
1702  IMPLICIT NONE
1703 
1704 #include "Indexing.h"
1705 
1706 ! ... parameters
1707  TYPE(t_region) :: region
1708  TYPE(t_patch) :: patch
1709 
1710 ! ... loop variables
1711  INTEGER :: idum, i, j, k
1712 
1713 ! ... local variables
1714  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, idir, jdir, kdir
1715  INTEGER :: ilev, lbound, icoff, ijcoff, ijkc, ijkd
1716  INTEGER :: gasmodel
1717 
1718  REAL(RFREAL) :: sgn(3)
1719  REAL(RFREAL), POINTER :: cv(:,:)
1720 
1721 !******************************************************************************
1722 
1723  CALL registerfunction( region%global,'RFLO_BcondSymmetry',&
1724  'RFLO_ModBoundaryConditions.F90' )
1725 
1726 ! get dimensions and pointers
1727 
1728  ilev = region%currLevel
1729  lbound = patch%lbound
1730 
1731  CALL rflo_getpatchindices( region,patch,ilev, &
1732  ibeg,iend,jbeg,jend,kbeg,kend )
1734  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
1735 
1736  gasmodel = region%mixtInput%gasModel
1737 
1738  cv => region%levels(ilev)%mixt%cv
1739 
1740 ! to mirror the appropriate velocity components
1741 
1742  IF (lbound==1 .OR. lbound==2) THEN
1743  sgn(1) = -1._rfreal
1744  sgn(2) = +1._rfreal
1745  sgn(3) = +1._rfreal
1746  ELSE IF (lbound==3 .OR. lbound==4) THEN
1747  sgn(1) = +1._rfreal
1748  sgn(2) = -1._rfreal
1749  sgn(3) = +1._rfreal
1750  ELSE IF (lbound==5 .OR. lbound==6) THEN
1751  sgn(1) = +1._rfreal
1752  sgn(2) = +1._rfreal
1753  sgn(3) = -1._rfreal
1754  ENDIF
1755 
1756 ! loop over all cells of a patch
1757 
1758  DO idum=1,region%nDumCells
1759  DO k=kbeg,kend
1760  DO j=jbeg,jend
1761  DO i=ibeg,iend
1762  ijkc = indijk(i+(idum-1)*idir,j+(idum-1)*jdir,k+(idum-1)*kdir,icoff,ijcoff)
1763  ijkd = indijk(i-idum*idir,j-idum*jdir,k-idum*kdir,icoff,ijcoff)
1764 
1765  cv(cv_mixt_dens,ijkd) = cv(cv_mixt_dens,ijkc)
1766  cv(cv_mixt_xmom,ijkd) = sgn(1)*cv(cv_mixt_xmom,ijkc)
1767  cv(cv_mixt_ymom,ijkd) = sgn(2)*cv(cv_mixt_ymom,ijkc)
1768  cv(cv_mixt_zmom,ijkd) = sgn(3)*cv(cv_mixt_zmom,ijkc)
1769  cv(cv_mixt_ener,ijkd) = cv(cv_mixt_ener,ijkc)
1770 
1771  IF (gasmodel == gas_model_tcperf) THEN
1772  CALL mixtureproperties( region,ijkd,ijkd,.false. )
1773  ELSE
1774  CALL mixtureproperties( region,ijkd,ijkd,.true. )
1775  ENDIF
1776  ENDDO ! i
1777  ENDDO ! j
1778  ENDDO ! k
1779  ENDDO ! idum
1780 
1781 ! finalize
1782 
1783  CALL deregisterfunction( region%global )
1784 
1785 END SUBROUTINE rflo_bcondsymmetry
1786 
1787 !******************************************************************************
1788 !
1789 ! Purpose: receive data for dummy cells from adjacent regions being
1790 ! on another processor.
1791 !
1792 ! Description: none.
1793 !
1794 ! Input: regions = data of all regions
1795 ! iReg = index of current region.
1796 !
1797 ! Output: regions(iReg)%levels%mixt = updated flow values (cv,dv,tv,gv)
1798 ! in dummy cells of current region.
1799 !
1800 ! Notes: none.
1801 !
1802 !******************************************************************************
1803 
1804 SUBROUTINE rflo_boundaryconditionsrecv( regions,iReg )
1805 
1809  IMPLICIT NONE
1810 
1811 ! ... parameters
1812  TYPE(t_region), POINTER :: regions(:)
1813 
1814  INTEGER :: ireg
1815 
1816 ! ... loop variables
1817  INTEGER :: ipatch
1818 
1819 ! ... local variables
1820  INTEGER :: ilev, npatches, bctype, iregsrc, ipatchsrc
1821 
1822  TYPE(t_patch), POINTER :: patch, patchsrc
1823 
1824 !******************************************************************************
1825 
1826  CALL registerfunction( regions(ireg)%global,'RFLO_BoundaryConditionsRecv', &
1827  'RFLO_ModBoundaryConditions.F90' )
1828 
1829 ! get dimensions --------------------------------------------------------------
1830 
1831  ilev = regions(ireg)%currLevel
1832  npatches = regions(ireg)%nPatches
1833 
1834 ! receive data (regular cells) from other processors --------------------------
1835 
1836  DO ipatch=1,npatches
1837 
1838  patch => regions(ireg)%levels(ilev)%patches(ipatch)
1839 
1840  bctype = patch%bcType
1841  iregsrc = patch%srcRegion
1842  ipatchsrc = patch%srcPatch
1843 
1844 ! - region interface, periodic boundary
1845 
1846  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
1847  (bctype>=bc_regionint .AND. bctype<=bc_regionint +bc_range) .OR. &
1848  (bctype>=bc_regnonconf .AND. bctype<=bc_regnonconf+bc_range) .OR. &
1849  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
1850  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
1851  IF (regions(iregsrc)%procid /= regions(ireg)%global%myProcid) THEN
1852  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1853 
1854  CALL rflo_receivedummyvals( regions(ireg),regions(iregsrc), &
1855  patch,patchsrc )
1856  ENDIF
1857  ENDIF
1858 
1859  ENDDO ! iPatch
1860 
1861 ! copy/exchange data for edge and corner cells --------------------------------
1862 
1863  CALL rflo_setcorneredgecells( regions(ireg) )
1864 
1865  DO ipatch=1,npatches
1866  patch => regions(ireg)%levels(ilev)%patches(ipatch)
1867  bctype = patch%bcType
1868  IF ((bctype>=bc_noslipwall .AND. bctype<=bc_noslipwall+bc_range) .OR. &
1869  (bctype>=bc_injection .AND. bctype<=bc_injection +bc_range .AND. &
1870  regions(ireg)%mixtInput%flowModel==flow_navst) .OR. &
1871  (bctype>=bc_symmetry .AND. bctype<=bc_symmetry +bc_range)) THEN
1872  CALL rflo_correctcorneredgecells( regions(ireg),patch,bctype )
1873  ENDIF
1874  ENDDO ! iPatch
1875 
1876  CALL rflo_exchangecorneredgecells( regions,ireg )
1877  CALL rflo_receivecorneredgecells( regions,ireg )
1878 
1879 ! finalize --------------------------------------------------------------------
1880 
1881  CALL deregisterfunction( regions(ireg)%global )
1882 
1883 END SUBROUTINE rflo_boundaryconditionsrecv
1884 
1885 !******************************************************************************
1886 !
1887 ! Purpose: send data to dummy cells of adjacent regions being located
1888 ! on another processor.
1889 !
1890 ! Description: none.
1891 !
1892 ! Input: regions = data of all regions
1893 ! iReg = index of current region.
1894 !
1895 ! Output: data to other processors.
1896 !
1897 ! Notes: none.
1898 !
1899 !******************************************************************************
1900 
1901 SUBROUTINE rflo_boundaryconditionssend( regions,iReg )
1902 
1905  IMPLICIT NONE
1906 
1907 ! ... parameters
1908  TYPE(t_region), POINTER :: regions(:)
1909 
1910  INTEGER :: ireg
1911 
1912 ! ... loop variables
1913  INTEGER :: ipatch
1914 
1915 ! ... local variables
1916  INTEGER :: ilev, npatches, bctype, iregsrc, ipatchsrc
1917 
1918  TYPE(t_patch), POINTER :: patch, patchsrc
1919  TYPE(t_global), POINTER :: global
1920 
1921 !******************************************************************************
1922 
1923  global => regions(ireg)%global
1924 
1925  CALL registerfunction( global,'RFLO_BoundaryConditionsSend',&
1926  'RFLO_ModBoundaryConditions.F90' )
1927 
1928 ! get dimensions --------------------------------------------------------------
1929 
1930  ilev = regions(ireg)%currLevel
1931  npatches = regions(ireg)%nPatches
1932 
1933 ! send data for edge and corner cells -----------------------------------------
1934 
1935  CALL rflo_sendcorneredgecells( regions,ireg )
1936 
1937 ! loop over patches -----------------------------------------------------------
1938 
1939  DO ipatch=1,npatches
1940 
1941  patch => regions(ireg)%levels(ilev)%patches(ipatch)
1942 
1943  bctype = patch%bcType
1944  iregsrc = patch%srcRegion
1945  ipatchsrc = patch%srcPatch
1946 
1947 ! - conforming region interface
1948 
1949  IF (bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) THEN
1950  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1951 
1952  IF (regions(iregsrc)%procid /= global%myProcid) THEN
1953  CALL rflo_senddummyconf( regions(ireg),regions(iregsrc),patch )
1954  ENDIF
1955 
1956 ! - non-conforming region interface (integer)
1957 
1958  ELSE IF (bctype>=bc_regionint .AND. bctype<=bc_regionint+bc_range) THEN
1959  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1960 
1961  IF (regions(iregsrc)%procid /= global%myProcid) THEN
1962  CALL rflo_senddummyint( regions(ireg),regions(iregsrc),patch )
1963  ENDIF
1964 
1965 ! - non-conforming region interface (irregular)
1966 
1967  ELSE IF (bctype>=bc_regnonconf .AND. bctype<=bc_regnonconf+bc_range) THEN
1968  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1969 
1970  IF (regions(iregsrc)%procid /= global%myProcid) THEN
1971  CALL rflo_senddummyireg( regions(ireg),regions(iregsrc),patch )
1972  ENDIF
1973 
1974 ! - translational periodicity
1975 
1976  ELSE IF (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri+bc_range) THEN
1977  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1978 
1979  IF (regions(iregsrc)%procid /= global%myProcid) THEN
1980  CALL rflo_senddummyconf( regions(ireg),regions(iregsrc),patch )
1981  ENDIF
1982 
1983 ! - rotational periodicity
1984 
1985  ELSE IF (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri+bc_range) THEN
1986  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
1987 
1988  CALL rflo_bcondrotatperiod( regions(ireg),regions(iregsrc), &
1989  patch,patchsrc )
1990 
1991  ENDIF ! bcType
1992 
1993  ENDDO ! iPatch
1994 
1995 ! finalize --------------------------------------------------------------------
1996 
1997  CALL deregisterfunction( global )
1998 
1999 END SUBROUTINE rflo_boundaryconditionssend
2000 
2001 !******************************************************************************
2002 !
2003 ! Purpose: set boundary conditions or exchange data between adjacent
2004 ! regions being on same processor.
2005 !
2006 ! Description: none.
2007 !
2008 ! Input: regions = data of all regions
2009 ! iReg = index of current region.
2010 !
2011 ! Output: regions(iReg)%levels%mixt = updated flow values (cv,dv,tv,gv)
2012 ! in dummy cells of current region.
2013 !
2014 ! Notes: none.
2015 !
2016 !******************************************************************************
2017 
2018 SUBROUTINE rflo_boundaryconditionsset( regions,iReg )
2019 
2022  IMPLICIT NONE
2023 
2024 ! ... parameters
2025  TYPE(t_region), POINTER :: regions(:)
2026 
2027  INTEGER :: ireg
2028 
2029 ! ... loop variables
2030  INTEGER :: ipatch
2031 
2032 ! ... local variables
2033  INTEGER :: ilev, npatches, bctype, iregsrc, ipatchsrc
2034 
2035  TYPE(t_patch), POINTER :: patch, patchsrc
2036  TYPE(t_global), POINTER :: global
2037 
2038 !******************************************************************************
2039 
2040  global => regions(ireg)%global
2041 
2042  CALL registerfunction( global,'RFLO_BoundaryConditionsSet',&
2043  'RFLO_ModBoundaryConditions.F90' )
2044 
2045 ! get dimensions --------------------------------------------------------------
2046 
2047  ilev = regions(ireg)%currLevel
2048  npatches = regions(ireg)%nPatches
2049 
2050 ! loop over patches -----------------------------------------------------------
2051 
2052  DO ipatch=1,npatches
2053 
2054  patch => regions(ireg)%levels(ilev)%patches(ipatch)
2055 
2056  bctype = patch%bcType
2057  iregsrc = patch%srcRegion
2058  ipatchsrc = patch%srcPatch
2059 
2060 ! - inflow
2061 
2062  IF (bctype>=bc_inflow .AND. bctype<=bc_inflow+bc_range) THEN
2063  CALL rflo_bcondinflow( regions(ireg),patch )
2064 
2065 ! - outflow
2066 
2067  ELSE IF (bctype>=bc_outflow .AND. bctype<=bc_outflow+bc_range) THEN
2068  CALL rflo_bcondoutflow( regions(ireg),patch )
2069 
2070 ! - conforming region interface
2071 
2072  ELSE IF (bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) THEN
2073  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
2074 
2075  IF (regions(iregsrc)%procid == global%myProcid) THEN
2076  CALL rflo_exchangedummyconf( regions(ireg),regions(iregsrc), &
2077  patch,patchsrc )
2078  ENDIF
2079 
2080 ! - non-conforming region interface (integer)
2081 
2082  ELSE IF (bctype>=bc_regionint .AND. bctype<=bc_regionint+bc_range) THEN
2083  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
2084 
2085  IF (regions(iregsrc)%procid == global%myProcid) THEN
2086  CALL rflo_exchangedummyint( regions(ireg),regions(iregsrc), &
2087  patch,patchsrc )
2088  ENDIF
2089 
2090 ! - non-conforming region interface (irregular)
2091 
2092  ELSE IF (bctype>=bc_regnonconf .AND. bctype<=bc_regnonconf+bc_range) THEN
2093  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
2094 
2095  IF (regions(iregsrc)%procid == global%myProcid) THEN
2096  CALL rflo_exchangedummyireg( regions(ireg),regions(iregsrc), &
2097  patch,patchsrc )
2098  ENDIF
2099 
2100 ! - slip wall
2101 
2102  ELSE IF (bctype>=bc_slipwall .AND. bctype<=bc_slipwall+bc_range) THEN
2103  CALL rflo_bcondslipwall( regions(ireg),patch )
2104 
2105 ! - noslip wall
2106 
2107  ELSE IF (bctype>=bc_noslipwall .AND. bctype<=bc_noslipwall+bc_range) THEN
2108  CALL rflo_bcondnoslipwall( regions(ireg),patch )
2109 
2110 ! - far field
2111 
2112  ELSE IF (bctype>=bc_farfield .AND. bctype<=bc_farfield+bc_range) THEN
2113  CALL rflo_bcondfarfield( regions(ireg),patch )
2114 
2115 ! - injection
2116 
2117  ELSE IF (bctype>=bc_injection .AND. bctype<=bc_injection+bc_range) THEN
2118  IF (bctype==bc_injection_apn) THEN
2119  IF (patch%bcCoupled==bc_external) THEN
2120  CALL errorstop( global,err_external_funct,&
2121  __line__, &
2122  'BC_INJECT_APN used in .bc with interaction flag in .top file' )
2123  ELSE
2124  CALL rflo_bcondinjectionapn( regions(ireg),patch )
2125  ENDIF
2126  ENDIF
2127  CALL rflo_bcondinjection( regions(ireg),patch )
2128 
2129 ! - symmetry
2130 
2131  ELSE IF (bctype>=bc_symmetry .AND. bctype<=bc_symmetry+bc_range) THEN
2132  CALL rflo_bcondsymmetry( regions(ireg),patch )
2133 
2134 ! - translational periodicity
2135 
2136  ELSE IF (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri+bc_range) THEN
2137  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
2138 
2139  IF (regions(iregsrc)%procid == global%myProcid) THEN
2140  CALL rflo_exchangedummyconf( regions(ireg),regions(iregsrc), &
2141  patch,patchsrc )
2142  ENDIF
2143 
2144 ! - rotational periodicity
2145 
2146  ELSE IF (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri+bc_range) THEN
2147  patchsrc => regions(iregsrc)%levels(ilev)%patches(ipatchsrc)
2148 
2149  CALL rflo_bcondrotatperiod( regions(ireg),regions(iregsrc), &
2150  patch,patchsrc )
2151 
2152 ! - I dunno ...
2153 
2154  ELSE
2155  CALL errorstop( global,err_unknown_bc,&
2156  __line__ )
2157  ENDIF ! bcType
2158 
2159  ENDDO ! iPatch
2160 
2161 ! finalize --------------------------------------------------------------------
2162 
2163  CALL deregisterfunction( global )
2164 
2165 END SUBROUTINE rflo_boundaryconditionsset
2166 
2167 !******************************************************************************
2168 !
2169 ! Purpose: set outflow boundary condition for one cell.
2170 !
2171 ! Description: the subsonic boundary condition is based on non-reflecting,
2172 ! method of Rudy and Strikwerda: Boundary Conditions for Subsonic
2173 ! Compressible Navier-Stokes Calculations.
2174 ! Computers and Fluids 9, 327-338 (1981). The supersonic boundary
2175 ! condition consists of simple extrapolation.
2176 !
2177 ! Input: bcOpt = boundary treatment: subsonic, supersonic, or mixed
2178 ! pout = given static outlet pressure
2179 ! sx/y/zn = components of ortho-normalized face vector (outward facing)
2180 ! cpgas = specific heat at constant pressure (boundary cell)
2181 ! mol = molecular mass at boundary cell
2182 ! rho = density at boundary cell
2183 ! rhou/v/w = density * velocity components at boundary cell
2184 ! rhoe = density * total energy at boundary cell
2185 ! press = static pressure at boundary cell
2186 !
2187 ! Output: rhob = density at boundary
2188 ! rhou/v/wb = density * velocity components at boundary
2189 ! rhoeb = density * total energy at boundary
2190 !
2191 ! Notes: this condition is valid only for thermally and calorically
2192 ! perfect gas (supersonic outflow valid for all gases).
2193 !
2194 !******************************************************************************
2195 !
2196 ! $Id: RFLO_ModBoundaryConditions.F90,v 1.30 2010/02/18 21:47:39 juzhang Exp $
2197 !
2198 ! Copyright: (c) 2002 by the University of Illinois
2199 !
2200 !******************************************************************************
2201 
2202 SUBROUTINE rflo_bcondoutflowperf( bcOpt,pout,sxn,syn,szn,cpgas,mol, &
2203  rho,rhou,rhov,rhow,rhoe,press, &
2204  ud,vd,wd,tempd,uold,vold,wold, udold,vdold, &
2205  wdold,rdold,edold, reflength, dtmin, kappa, &
2206  rhob,rhoub,rhovb,rhowb,rhoeb )
2207 
2208  USE moddatatypes
2209  USE modparameters
2212 
2213  IMPLICIT NONE
2214 
2215 ! ... parameters
2216  INTEGER :: bcopt
2217 
2218  REAL(RFREAL) :: pout
2219  REAL(RFREAL) :: rho, rhou, rhov, rhow, rhoe, press
2220  REAL(RFREAL) :: ud, vd, wd, tempd, uold, vold, wold
2221  REAL(RFREAL) :: udold, vdold, wdold, rdold, edold, reflength, dtmin, kappa
2222  REAL(RFREAL) :: sxn, syn, szn, csound, cpgas, mol
2223  REAL(RFREAL) :: rhob, rhoub, rhovb, rhowb, rhoeb
2224 
2225 ! ... local variables
2226  REAL(RFREAL) :: rgas, gamma, gam1, u, v, w, mach, rhoc, deltp, &
2227  ub, vb, wb, vnd
2228  REAL(RFREAL) :: pdold, qdold2, qnold, qn, pb, beta
2229 
2230 !******************************************************************************
2231 ! gas properties; velocity components; Mach number
2232 
2233  rgas = mixtperf_r_m( mol )
2234  gamma = mixtperf_g_cpr( cpgas,rgas )
2235  gam1 = gamma - 1.0_rfreal
2236 
2237  u = rhou/rho
2238  v = rhov/rho
2239  w = rhow/rho
2240  csound = mixtperf_c_dgp( rho,gamma,press )
2241  mach = sqrt(u*u+v*v+w*w)/csound
2242 
2243 ! subsonic flow ---------------------------------------------------------------
2244 
2245  IF (mach < 1.0_rfreal .AND. &
2246  (bcopt == bcopt_subsonic .OR. bcopt == bcopt_mixed)) THEN
2247 
2248  beta = kappa*csound*(1._rfreal - mach*mach)/reflength
2249  rhoc = rho*csound
2250  qdold2= udold*udold + vdold*vdold + wdold*wdold
2251  pdold = mixtperf_p_deogvm2( rdold,edold,gamma,qdold2 )
2252  deltp = pdold - pout
2253  qn = u*sxn + v*syn + w*szn
2254  qnold = uold*sxn + vold*syn + wold*szn
2255 ! qn = ud*sxn + vd*syn + wd*szn
2256 ! qnOld = udOld*sxn + vdOld*syn + wdOld*szn
2257 ! qn = 0.5_RFREAL*((ud+u)*sxn + (vd+v)*syn + (wd+w)*szn)
2258 ! qnOld = 0.5_RFREAL*((udOld+uOld)*sxn + (vdOld+vOld)*syn + (wdOld+wOld)*szn)
2259 
2260  pb = press + rhoc*(qn - qnold) - dtmin*beta*(press - pout)
2261 ! pb = pdOld + rhoc*(qn - qnOld) - dtMin*beta*(pdOld - pout)
2262 ! pb = pout
2263  ub = ud
2264  vb = vd
2265  wb = wd
2266  rhob = pb/(rgas*tempd)
2267 
2268 ! - special treatment to prevent "deltp" from changing the sign
2269 ! of velocity components. This may happen for very small u, v, w.
2270 
2271  vnd = ub*sxn + vb*syn + wb*szn
2272  IF ( vnd < 0.0_rfreal ) THEN ! inflow at outflow boundary
2273  ub = 0.010_rfreal*sign(1.0_rfreal,u)*min(abs(ub),abs(u))
2274  vb = 0.010_rfreal*sign(1.0_rfreal,v)*min(abs(vb),abs(v))
2275  wb = 0.010_rfreal*sign(1.0_rfreal,w)*min(abs(wb),abs(w))
2276  END IF ! vnd
2277 
2278  rhoub = rhob*ub
2279  rhovb = rhob*vb
2280  rhowb = rhob*wb
2281  rhoeb = rhob*mixtperf_eo_dgpuvw( rhob,gamma,pb,ub,vb,wb )
2282 
2283 ! supersonic flow -------------------------------------------------------------
2284 
2285  ELSE
2286  rhob = rho
2287  rhoub = rhou
2288  rhovb = rhov
2289  rhowb = rhow
2290  rhoeb = rhoe
2291  END IF ! mach
2292 
2293 END SUBROUTINE rflo_bcondoutflowperf
2294 
2295 !******************************************************************************
2296 !
2297 ! Purpose: set inflow boundary condition for one cell based on presribed
2298 ! velocity and temperature.
2299 !
2300 ! Description: the subsonic boundary condition is based on non-reflecting,
2301 ! method of Rudy and Strikwerda: Boundary Conditions for Subsonic
2302 ! Compressible Navier-Stokes Calculations.
2303 ! Computers and Fluids 9, 327-338 (1981). The supersonic boundary
2304 ! condition consists of simple extrapolation.
2305 !
2306 ! Input: bcOpt = inflow bc-option: steady or unsteady
2307 ! uInfl = prescribed u-velocity
2308 ! vInfl = prescribed v-velocity
2309 ! wInfl = prescribed w-velocity
2310 ! tInfl = prescribed temperature
2311 ! pInfl = prescribed pressure (only used if supersonic)
2312 ! sx/y/zn = components of ortho-normalized face vector (outward facing)
2313 ! cpgas = specific heat at constant pressure (boundary cell)
2314 ! mol = molecular mass at boundary cell
2315 ! rho = density at boundary cell
2316 ! rhou/v/w = density * velocity components at boundary cell
2317 ! rhoe = density * total energy at boundary cell
2318 ! press = static pressure at boundary cell
2319 !
2320 ! Output: rhob = density at boundary
2321 ! rhou/v/wb = density * velocity components at boundary
2322 ! rhoeb = density * total energy at boundary
2323 !
2324 ! Notes: this condition is valid only for a thermally and calorically
2325 ! perfect gas (supersonic in/outflow valid for all gases).
2326 !
2327 !******************************************************************************
2328 
2329 SUBROUTINE rflo_bcondinflowvtperf( model,bcOpt,uInfl,vInfl,wInfl,tInfl,pInfl, &
2330  sxn,syn,szn,cpgas,mol, &
2331  rho,rhou,rhov,rhow,rhoe,press, &
2332  rhob,rhoub,rhovb,rhowb,rhoeb )
2333 
2336 
2337  IMPLICIT NONE
2338 
2339 ! ... parameters
2340  INTEGER :: model, bcopt
2341  REAL(RFREAL) :: uinfl, vinfl, winfl, tinfl, pinfl
2342  REAL(RFREAL) :: sxn, syn, szn, cpgas, mol
2343  REAL(RFREAL) :: rho, rhou, rhov, rhow, rhoe, press
2344  REAL(RFREAL) :: rhob, rhoub, rhovb, rhowb, rhoeb
2345 
2346 ! ... local variables
2347  INTEGER :: ierror
2348  REAL(RFREAL) :: rgas, gamma, mach
2349  REAL(RFREAL) :: re, ue, ve, we, pe, qn, crho0, csound
2350  REAL(RFREAL) :: ub, vb, wb, tb, pb
2351 
2352 !******************************************************************************
2353 ! gas properties
2354 
2355  rgas = mixtperf_r_m( mol )
2356  gamma = mixtperf_g_cpr( cpgas,rgas )
2357 
2358 ! flow values at a reference location (= interior cell)
2359 
2360  re = rho
2361  ue = rhou/rho
2362  ve = rhov/rho
2363  we = rhow/rho
2364  pe = press
2365 
2366 ! flow values at inlet plane
2367 
2368  tb = tinfl
2369  ub = uinfl
2370  vb = vinfl
2371  wb = winfl
2372  qn = sxn*ub + syn*vb + szn*wb
2373  mach= qn/mixtperf_c_grt( gamma,rgas,tb )
2374 
2375  IF (qn > 0 .AND. model==bcopt_steady) THEN
2376  WRITE(stderr,*)'BC_INFLOW_VELTEMP: outflow detected at inflow plane'
2377 #ifdef MPI
2378  CALL mpi_abort( ierror )
2379 #endif
2380  stop
2381  ENDIF
2382 
2383  IF (bcopt == bcopt_mixed) THEN
2384 
2385 ! - subsonic inflow (qn<0)
2386 
2387  IF (mach < 1._rfreal) THEN
2388  csound = mixtperf_c_dgp( re,gamma,pe )
2389  crho0 = csound*re
2390 
2391  pb = pe-crho0*(sxn*(ub-ue)+syn*(vb-ve)+szn*(wb-we))
2392 
2393 ! - supersonic inflow (qn<0)
2394 
2395  ELSE
2396 
2397  pb = pinfl
2398  ENDIF
2399 
2400  ELSEIF (bcopt == bcopt_subsonic) THEN
2401  csound = mixtperf_c_dgp( re,gamma,pe )
2402  crho0 = csound*re
2403 
2404  pb = pe-crho0*(sxn*(ub-ue)+syn*(vb-ve)+szn*(wb-we))
2405 
2406  ELSEIF (bcopt == bcopt_supersonic) THEN
2407 
2408  pb = pinfl
2409  ENDIF
2410 
2411  rhob = mixtperf_d_prt( pb,rgas,tb )
2412  rhoub = rhob*ub
2413  rhovb = rhob*vb
2414  rhowb = rhob*wb
2415  rhoeb = rhob*mixtperf_eo_dgpuvw( rhob,gamma,pb,ub,vb,wb )
2416 
2417 END SUBROUTINE rflo_bcondinflowvtperf
2418 
2419 !******************************************************************************
2420 !
2421 ! Purpose: set inflow boundary condition for one cell based on presribed
2422 ! velocity and pressure.
2423 !
2424 ! Description: the subsonic boundary condition is based on non-reflecting,
2425 ! method of Rudy and Strikwerda: Boundary Conditions for Subsonic
2426 ! Compressible Navier-Stokes Calculations.
2427 ! Computers and Fluids 9, 327-338 (1981). The supersonic boundary
2428 ! condition consists of simple extrapolation.
2429 !
2430 ! Input: bcOpt = inflow bc-option: steady or unsteady
2431 ! uInfl = prescribed u-velocity
2432 ! vInfl = prescribed v-velocity
2433 ! wInfl = prescribed w-velocity
2434 ! pInfl = prescribed pressure
2435 ! tInfl = prescribed temperature (only used if supersonic)
2436 ! sx/y/zn = components of ortho-normalized face vector (outward facing)
2437 ! cpgas = specific heat at constant pressure (boundary cell)
2438 ! mol = molecular mass at boundary cell
2439 ! rho = density at boundary cell
2440 ! rhou/v/w = density * velocity components at boundary cell
2441 ! rhoe = density * total energy at boundary cell
2442 ! press = static pressure at boundary cell
2443 !
2444 ! Output: rhob = density at boundary
2445 ! rhou/v/wb = density * velocity components at boundary
2446 ! rhoeb = density * total energy at boundary
2447 !
2448 ! Notes: this condition is valid only for a thermally and calorically
2449 ! perfect gas (supersonic in/outflow valid for all gases).
2450 !
2451 !******************************************************************************
2452 
2453 SUBROUTINE rflo_bcondinflowvpperf( model,bcOpt,uInfl,vInfl,wInfl,tInfl,pInfl, &
2454  sxn,syn,szn,cpgas,mol, &
2455  rho,rhou,rhov,rhow,rhoe,press, &
2456  rhob,rhoub,rhovb,rhowb,rhoeb )
2457 
2460 
2461  IMPLICIT NONE
2462 
2463 ! ... parameters
2464  INTEGER :: model, bcopt
2465  REAL(RFREAL) :: uinfl, vinfl, winfl, tinfl, pinfl
2466  REAL(RFREAL) :: sxn, syn, szn, cpgas, mol
2467  REAL(RFREAL) :: rho, rhou, rhov, rhow, rhoe, press
2468  REAL(RFREAL) :: rhob, rhoub, rhovb, rhowb, rhoeb
2469 
2470 ! ... local variables
2471  INTEGER :: ierror
2472  REAL(RFREAL) :: rgas, gamma, mach
2473  REAL(RFREAL) :: re, ue, ve, we, pe, te, qn, crho0, csound
2474  REAL(RFREAL) :: ub, vb, wb, tb, pb, ta
2475 
2476 !******************************************************************************
2477 ! gas properties
2478 
2479  rgas = mixtperf_r_m( mol )
2480  gamma = mixtperf_g_cpr( cpgas,rgas )
2481 
2482 ! flow values at a reference location (= interior cell)
2483 
2484  re = rho
2485  ue = rhou/rho
2486  ve = rhov/rho
2487  we = rhow/rho
2488  pe = press
2489  te = pe/(re*rgas)
2490 
2491 ! flow values at inlet plane
2492 
2493  pb = pinfl
2494  ub = uinfl
2495  vb = vinfl
2496  wb = winfl
2497  qn = sxn*ub + syn*vb + szn*wb
2498  ta = pb/(re*rgas)
2499  mach= qn/mixtperf_c_grt( gamma,rgas,ta )
2500 
2501  IF (qn > 0 .AND. model==bcopt_steady) THEN
2502  WRITE(stderr,*)'BC_INFLOW_VELPRESS: outflow detected at inflow plane'
2503 #ifdef MPI
2504  CALL mpi_abort( ierror )
2505 #endif
2506  stop
2507  ENDIF
2508 
2509  IF (bcopt == bcopt_mixed) THEN
2510 
2511 ! - subsonic inflow (qn<0)
2512 
2513  IF (mach < 1._rfreal) THEN
2514  csound = mixtperf_c_dgp( re,gamma,pe )
2515  crho0 = csound*re
2516 
2517  tb = te-csound/rgas*(sxn*(ub-ue)+syn*(vb-ve)+szn*(wb-we))
2518 
2519 ! - supersonic inflow (qn<0)
2520 
2521  ELSE
2522 
2523  tb = tinfl
2524  ENDIF
2525 
2526  ELSEIF (bcopt == bcopt_subsonic) THEN
2527  csound = mixtperf_c_dgp( re,gamma,pe )
2528  crho0 = csound*re
2529 
2530  tb = te-csound/rgas*(sxn*(ub-ue)+syn*(vb-ve)+szn*(wb-we))
2531 
2532  ELSEIF (bcopt == bcopt_supersonic) THEN
2533 
2534  tb = tinfl
2535  ENDIF
2536 
2537  rhob = mixtperf_d_prt( pb,rgas,tb )
2538  rhoub = rhob*ub
2539  rhovb = rhob*vb
2540  rhowb = rhob*wb
2541  rhoeb = rhob*mixtperf_eo_dgpuvw( rhob,gamma,pb,ub,vb,wb )
2542 
2543 END SUBROUTINE rflo_bcondinflowvpperf
2544 
2545 !******************************************************************************
2546 !
2547 ! Purpose: obtain fluctuations to be recycled at inflow boundary.
2548 !
2549 ! Description: take fluctuations = instantaneous - time-averaged-quantities,
2550 ! then add the fluctuations to corresp. imposed inlet quantities.
2551 !
2552 ! Input: bcOptModel = inflow-bc model: steady (laminar) or unsteady (turbulent)
2553 ! uvel = instantaneous u-velocity
2554 ! vvel = instantaneous v-velocity
2555 ! wvel = instantaneous w-velocity
2556 ! temp = instantaneous temperature (only used if supersonic)
2557 ! pres = instantaneous pressure
2558 ! muvel = time-averaged u-velocity
2559 ! mvvel = time-averaged v-velocity
2560 ! mwvel = time-averaged w-velocity
2561 ! mtemp = time-averaged temperature (only used if supersonic)
2562 ! mpres = time-averaged pressure
2563 !
2564 ! Output: ufluc = u fluctuation
2565 ! vfluc = v fluctuation
2566 ! wfluc = w fluctuation
2567 ! Tfluc = temperature fluctuation
2568 ! pfluc = pressure fluctuation
2569 !******************************************************************************
2570 SUBROUTINE rflo_bcondinflowfluc( bcOptModel, ufluc,vfluc,wfluc,Tfluc,pfluc, &
2571  uvel,vvel,wvel,temp,pres, &
2572  muvel,mvvel,mwvel,mtemp,mpres )
2573  IMPLICIT NONE
2574 
2575 ! ... parameter
2576  INTEGER :: bcoptmodel
2577  REAL(RFREAL) :: ufluc, vfluc, wfluc, tfluc, pfluc
2578  REAL(RFREAL) :: uvel, vvel, wvel, temp, pres
2579  REAL(RFREAL) :: muvel, mvvel, mwvel, mtemp, mpres
2580 
2581 ! ... local variables
2582 
2583 !******************************************************************************
2584 
2585 ! obtain fluctuations to be recycled
2586 
2587  IF (bcoptmodel==bcopt_unsteady) THEN
2588  ufluc = uvel-muvel
2589  vfluc = vvel-mvvel
2590  wfluc = wvel-mwvel
2591  tfluc = temp-mtemp
2592  pfluc = pres-mpres
2593  ELSE
2594  ufluc = 0._rfreal
2595  vfluc = 0._rfreal
2596  wfluc = 0._rfreal
2597  tfluc = 0._rfreal
2598  pfluc = 0._rfreal
2599  ENDIF
2600 
2601 END SUBROUTINE rflo_bcondinflowfluc
2602 
2603 ! ******************************************************************************
2604 ! End
2605 ! ******************************************************************************
2606 
2607 END MODULE rflo_modboundaryconditions
2608 
2609 ! ******************************************************************************
2610 !
2611 ! RCS Revision history:
2612 !
2613 ! $Log: RFLO_ModBoundaryConditions.F90,v $
2614 ! Revision 1.30 2010/02/18 21:47:39 juzhang
2615 ! Heat transfer bc for non-propellant surface documented in Rocburn_PY_HT.pdf in Rocburn_PY directory is implemented within Rocburn_PY. Major changes were made to Rocburn, Rocman3, RocfluidMP/genx, RocfluidMP/modflo directories.
2616 !
2617 ! Revision 1.29 2008/12/06 08:44:14 mtcampbe
2618 ! Updated license.
2619 !
2620 ! Revision 1.28 2008/11/19 22:17:27 mtcampbe
2621 ! Added Illinois Open Source License/Copyright
2622 !
2623 ! Revision 1.27 2006/08/19 15:38:50 mparmar
2624 ! Renamed patch variables
2625 !
2626 ! Revision 1.26 2006/05/05 02:19:32 wasistho
2627 ! commented redundant kernels in RFLO_BondInflow
2628 !
2629 ! Revision 1.25 2006/03/24 05:58:41 wasistho
2630 ! activated all components inlet fluctuations
2631 !
2632 ! Revision 1.24 2006/03/09 21:57:22 wasistho
2633 ! lower backflow coefficients in NR outflow bc
2634 !
2635 ! Revision 1.23 2006/02/18 08:28:10 wasistho
2636 ! added bcOpt sub/supersonic/mixed in veltemp and velpress
2637 !
2638 ! Revision 1.22 2006/01/30 20:16:20 wasistho
2639 ! added safety for using injectionAPN
2640 !
2641 ! Revision 1.21 2006/01/25 01:36:08 wasistho
2642 ! compute rhofVf in injectionAPN and mod injectionInit
2643 !
2644 ! Revision 1.20 2006/01/23 22:46:46 wasistho
2645 ! added bc_internal condition for injectionAPN
2646 !
2647 ! Revision 1.19 2006/01/21 11:38:10 wasistho
2648 ! added injectionapn
2649 !
2650 ! Revision 1.18 2005/12/07 02:25:15 wasistho
2651 ! added ifdef STATS in inflow veltemp and velpress
2652 !
2653 ! Revision 1.17 2005/12/06 22:08:22 wasistho
2654 ! ijkF to ijkD in the last fix
2655 !
2656 ! Revision 1.16 2005/12/06 21:54:03 wasistho
2657 ! control fluctuations for iflow veltemp and velpress
2658 !
2659 ! Revision 1.15 2005/12/01 09:33:43 wasistho
2660 ! added condition infloNijk < NIJK_INFLOW_INIT in bcondInflow
2661 !
2662 ! Revision 1.14 2005/11/30 23:36:06 wasistho
2663 ! increased tollerance in eps
2664 !
2665 ! Revision 1.13 2005/11/30 20:07:14 wasistho
2666 ! compute inflow fluctuations only when stats time>0
2667 !
2668 ! Revision 1.12 2005/11/29 22:47:52 wasistho
2669 ! bug fixed divided tav by integrTime in BcondInflow
2670 !
2671 ! Revision 1.11 2005/11/18 07:20:11 wasistho
2672 ! multiply veltemp/pres-inflow flucts by amplitude
2673 !
2674 ! Revision 1.10 2005/11/07 21:23:27 wasistho
2675 ! increased backflow prevention parameter
2676 !
2677 ! Revision 1.9 2005/10/31 22:34:04 wasistho
2678 ! added bcopt_steady in detecting inlet backflow
2679 !
2680 ! Revision 1.8 2005/10/31 21:09:35 haselbac
2681 ! Changed specModel and SPEC_MODEL_NONE
2682 !
2683 ! Revision 1.7 2005/10/26 22:23:31 wasistho
2684 ! lower backflow prevention parameter
2685 !
2686 ! Revision 1.6 2005/09/28 18:06:18 wasistho
2687 ! modified for moving slipwall
2688 !
2689 ! Revision 1.5 2005/09/20 23:13:29 wasistho
2690 ! added fluctuation option for inflow velTemp and velPress
2691 !
2692 ! Revision 1.4 2005/08/31 05:36:34 wasistho
2693 ! undefined nrcoef for supersonic
2694 !
2695 ! Revision 1.3 2005/05/13 08:12:38 wasistho
2696 ! implemented new bc_outflow model based on non-reflecting condition
2697 !
2698 ! Revision 1.2 2005/04/28 05:43:44 wasistho
2699 ! added velocity based inflow BC
2700 !
2701 ! Revision 1.1 2004/12/28 22:51:20 wasistho
2702 ! moved RFLO_Bcond* and RFLO_BoundaryCond* routines into RFLO_ModBoundaryConditions
2703 !
2704 !
2705 ! ******************************************************************************
2706 
2707 
2708 
2709 
2710 
2711 
2712 
2713 
2714 
2715 
2716 
2717 
2718 
2719 
2720 
2721 
2722 
2723 
2724 
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ibeg
subroutine rflo_bcondinjectionapn(region, patch)
real(rfreal) function mixtperf_p_deogvm2(D, Eo, G, Vm2)
Definition: MixtPerf_P.F90:39
subroutine rflo_getpatchdirection(patch, idir, jdir, kdir)
real(rfreal) function mixtperf_r_m(M)
Definition: MixtPerf_R.F90:54
subroutine bcondinjectionperf(distrib, minj, tinj, rhoVrel, sxn, syn, szn, cpgas, mm, p, rhob, rhoub, rhovb, rhowb, rhoeb, pb, uinj, vinj, winj)
static SURF_BEGIN_NAMESPACE double sign(double x)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE inode
subroutine bcondfarfieldperf(machInf, alphaInf, betaInf, pInf, tInf, sxn, syn, szn, cpgas, mol, rho, rhou, rhov, rhow, rhoe, press, rhob, rhoub, rhovb, rhowb, rhoeb, pb)
subroutine rflo_exchangedummyireg(region, regionSrc, patch, patchSrc)
subroutine, public rflo_bcondsymmetry(region, patch)
subroutine rflo_setcorneredgecells(region)
j indices k indices k
Definition: Indexing.h:6
subroutine rflo_sendcorneredgecells(regions, iReg)
subroutine, public rflo_bcondslipwall(region, patch)
real(rfreal) function mixtperf_c_dgp(D, G, P)
Definition: MixtPerf_C.F90:56
subroutine, public rflo_bcondinjection(region, patch)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine rflo_getpatchindices(region, patch, iLev, ibeg, iend, jbeg, jend, kbeg, kend)
subroutine rflo_exchangedummyconf(region, regionSrc, patch, patchSrc)
real(rfreal) function mixtperf_d_prt(P, R, T)
Definition: MixtPerf_D.F90:71
subroutine, public rflo_bcondrotatperiod(region, regionSrc, patch, patchSrc)
double sqrt(double d)
Definition: double.h:73
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
subroutine rflo_bcondinflowfluc(bcOptModel, ufluc, vfluc, wfluc, Tfluc, pfluc, uvel, vvel, wvel, temp, pres, muvel, mvvel, mwvel, mtemp, mpres)
subroutine, public rflo_bcondinflow(region, patch)
subroutine rflo_exchangecorneredgecells(regions, iReg)
subroutine rflo_getnodeoffset(region, iLev, iNodeOffset, ijNodeOffset)
subroutine rflo_receivecorneredgecells(regions, iReg)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jdir
subroutine rflo_receivedummyvals(region, regionSrc, patch, patchSrc)
subroutine rflo_senddummyint(region, regionSrc, patch)
subroutine bcondoutflowperf(bcOpt, pout, sxn, syn, szn, cpgas, mol, rho, rhou, rhov, rhow, rhoe, press, rhob, rhoub, rhovb, rhowb, rhoeb)
Definition: patch.h:74
subroutine rflo_correctcorneredgecells(region, patch, bcType)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode iend
subroutine rflo_bcondinflowvtperf(model, bcOpt, uInfl, vInfl, wInfl, tInfl, pInfl, sxn, syn, szn, cpgas, mol, rho, rhou, rhov, rhow, rhoe, press, rhob, rhoub, rhovb, rhowb, rhoeb)
subroutine, public rflo_boundaryconditionssend(regions, iReg)
subroutine, public rflo_bcondoutflow(region, patch)
subroutine, public rflo_boundaryconditionsrecv(regions, iReg)
subroutine, public rflo_bcondnoslipwall(region, patch)
blockLoc i
Definition: read.cpp:79
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE idir
subroutine rflo_getcelloffset(region, iLev, iCellOffset, ijCellOffset)
subroutine rflo_senddummyireg(region, regionSrc, patch)
subroutine bcondinflowperf(bcOptType, bcOptFixed, ptot, ttot, betah, betav, mach, sxn, syn, szn, cpgas, mm, rl, rul, rvl, rwl, rr, rur, rvr, rwr, rer, pr)
subroutine rflo_exchangedummyint(region, regionSrc, patch, patchSrc)
Vector_n min(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:346
subroutine rflo_senddummyconf(region, regionSrc, patch)
real(rfreal) function mixtperf_c_grt(G, R, T)
Definition: MixtPerf_C.F90:86
j indices j
Definition: Indexing.h:6
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jend
subroutine, public rflo_boundaryconditionsset(regions, iReg)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE kdir
real(rfreal) function mixtperf_eo_dgpuvw(D, G, P, U, V, W)
Definition: MixtPerf_E.F90:40
subroutine, public rflo_bcondfarfield(region, patch)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jbeg
subroutine rflo_copyboundarydata(global, patchPrev, patch)
subroutine mixtureproperties(region, inBeg, inEnd, gasUpdate)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode kbeg
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine rflo_bcondinflowvpperf(model, bcOpt, uInfl, vInfl, wInfl, tInfl, pInfl, sxn, syn, szn, cpgas, mol, rho, rhou, rhov, rhow, rhoe, press, rhob, rhoub, rhovb, rhowb, rhoeb)
real(rfreal) function mixtperf_g_cpr(Cp, R)
Definition: MixtPerf_G.F90:39
subroutine rflo_bcondoutflowperf(bcOpt, pout, sxn, syn, szn, cpgas, mol, rho, rhou, rhov, rhow, rhoe, press, ud, vd, wd, tempd, uOld, vOld, wOld, udOld, vdOld, wdOld, rdOld, edOld, refLength, dtMin, kappa, rhob, rhoub, rhovb, rhowb, rhoeb)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jnode
subroutine, public rflo_bcondinjectioninit(region)