Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_InitFlowHardCode.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: Initialize flow field in a region using hard-coded values.
26 !
27 ! Description: None.
28 !
29 ! Input:
30 ! pRegion Pointer to region
31 !
32 ! Output: None.
33 !
34 ! Notes:
35 ! 1. This routine assumes a perfect gas.
36 !
37 ! ******************************************************************************
38 !
39 ! $Id: RFLU_InitFlowHardCode.F90,v 1.38 2008/12/06 08:44:56 mtcampbe Exp $
40 !
41 ! Copyright: (c) 2003-2006 by the University of Illinois
42 !
43 ! ******************************************************************************
44 
45 SUBROUTINE rflu_initflowhardcode(pRegion)
46 
47  USE moddatatypes
48  USE moderror
49  USE moddatastruct, ONLY: t_region
50  USE modmixture, ONLY: t_mixt_input
51  USE modglobal, ONLY: t_global
52  USE modgrid, ONLY: t_grid
53  USE modparameters
54 
55  USE rflu_modbessel
58 
59  USE modinterfaces, ONLY: mixtperf_cv_cpr, &
66  mixtperf_r_cpg, &
68 
69  IMPLICIT NONE
70 
71 ! ******************************************************************************
72 ! Declarations and definitions
73 ! ******************************************************************************
74 
75 ! ==============================================================================
76 ! Arguments
77 ! ==============================================================================
78 
79  TYPE(t_region), POINTER :: pregion
80 
81 ! ==============================================================================
82 ! Locals
83 ! ==============================================================================
84 
85  CHARACTER(CHRLEN) :: errorstring,rcsidentstring
86  INTEGER :: ibc,icg,im,in,indcp,indmol,iq
87  REAL(RFREAL) :: atot,const,cp,cpref,cvm,cvg,cvl,cvv,d,dinc,dmin,doffs,dtot, &
88  dummyreal,etaqm,g,gc,gcref,gref,gx,gy,gz,h,height,l,mi,minj, &
89  mw,nx,ny,nz,omega,p,pg,pl,pmin,poffs,ptot,pv,r,rg,rgas,ri, &
90  rl,ro,rv,rvap,t,ttot,theta,u,ur,ut,uz,v,vinj,vm2,w,x,xmin, &
91  y,yp,ymin,z,zmin
92  REAL(RFREAL) :: xx,uo_c,a_c,l_l,c,a,po,uo,xo,rc_l,rc,a_cl,radius,psi,r1,r2,um
93  REAL(RFREAL), DIMENSION(:,:), POINTER :: pcv,pgv
94  TYPE(t_global), POINTER :: global
95  TYPE(t_grid), POINTER :: pgrid
96  TYPE(t_mixt_input), POINTER :: pmixtinput
97 
98 ! ******************************************************************************
99 ! Start
100 ! ******************************************************************************
101 
102  rcsidentstring = '$RCSfile: RFLU_InitFlowHardCode.F90,v $ $Revision: 1.38 $'
103 
104  global => pregion%global
105 
106  CALL registerfunction(global,'RFLU_InitFlowHardCode', &
107  'RFLU_InitFlowHardCode.F90')
108 
109  IF ( global%verbLevel > verbose_none ) THEN
110  WRITE(stdout,'(A,1X,A)') solver_name, &
111  'Initializing flow field from hard code...'
112 
113  IF ( global%verbLevel > verbose_low ) THEN
114  WRITE(stdout,'(A,3X,A,A)') solver_name,'Case: ',trim(global%casename)
115  END IF ! global%verbLevel
116  END IF ! global%verbLevel
117 
118 ! ******************************************************************************
119 ! Set pointers
120 ! ******************************************************************************
121 
122  pgrid => pregion%grid
123  pcv => pregion%mixt%cv
124  pgv => pregion%mixt%gv
125  pmixtinput => pregion%mixtInput
126 
127  indcp = pregion%mixtInput%indCp
128  indmol = pregion%mixtInput%indMol
129 
130 ! ******************************************************************************
131 ! Set constant gas properties. NOTE used only in some cases
132 ! ******************************************************************************
133 
134  cpref = global%refCp
135  gref = global%refGamma
136  gcref = mixtperf_r_cpg(cpref,gref)
137 
138 ! ******************************************************************************
139 ! Initialize flow field based on user input and fluid model
140 ! ******************************************************************************
141 
142  SELECT CASE ( pmixtinput%fluidModel )
143 
144 ! ==============================================================================
145 ! Incompressible fluid model
146 ! ==============================================================================
147 
148  CASE ( fluid_model_incomp )
149  pregion%mixt%cvState = cv_mixt_state_prim
150 
151 ! TEMPORARY, to be replaced by proper code
152  CALL errorstop(global,err_reached_default,__line__)
153 ! END TEMPORARY
154 
155 ! ==============================================================================
156 ! Compressible fluid model
157 ! ==============================================================================
158 
159  CASE ( fluid_model_comp )
160  pregion%mixt%cvState = cv_mixt_state_cons
161 
162  SELECT CASE ( global%casename )
163 
164 ! ------------------------------------------------------------------------------
165 ! Cylinder shock diffraction
166 ! ------------------------------------------------------------------------------
167 
168  CASE ( "cylds" )
169  DO icg = 1,pgrid%nCellsTot
170  x = pgrid%cofg(xcoord,icg)
171 
172  IF ( x < pmixtinput%prepRealVal1 ) THEN
173  d = pmixtinput%prepRealVal2
174  u = pmixtinput%prepRealVal3
175  v = 0.0_rfreal
176  w = 0.0_rfreal
177  p = pmixtinput%prepRealVal4
178  ELSE
179  d = pmixtinput%prepRealVal5
180  u = 0.0_rfreal
181  v = 0.0_rfreal
182  w = 0.0_rfreal
183  p = pmixtinput%prepRealVal6
184  END IF ! x
185 
186  mw = pgv(gv_mixt_mol,indmol*icg)
187  cp = pgv(gv_mixt_cp ,indcp *icg)
188 
189  gc = mixtperf_r_m(mw)
190  g = mixtperf_g_cpr(cp,gc)
191 
192  pcv(cv_mixt_dens,icg) = d
193  pcv(cv_mixt_xmom,icg) = d*u
194  pcv(cv_mixt_ymom,icg) = d*v
195  pcv(cv_mixt_zmom,icg) = d*w
196  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
197  END DO ! icg
198 
199 ! ------------------------------------------------------------------------------
200 ! Diffracting shock issuing from open-ended shock tube
201 ! ------------------------------------------------------------------------------
202 
203 ! ----- Without backchamber ----------------------------------------------------
204 
205  CASE ( "ds_7p25" , "ds_20p0" , "ds_50p0" , "ds_125p0" , &
206  "ds_7p25_v2", "ds_20p0_v2", "ds_50p0_v2", "ds_125p0_v2", &
207  "ds_7p25_v3", "ds_20p0_v3", "ds_50p0_v3", "ds_125p0_v3", &
208  "ds_7p25_v4", "ds_20p0_v4", "ds_50p0_v4", "ds_125p0_v4", &
209  "ds_7p25_v5", "ds_20p0_v5", "ds_50p0_v5", "ds_125p0_v5" )
210  DO icg = 1,pgrid%nCellsTot
211  x = pgrid%cofg(xcoord,icg)
212 
213  IF ( x < pmixtinput%prepRealVal1 ) THEN
214  d = pmixtinput%prepRealVal2
215  u = pmixtinput%prepRealVal3
216  v = 0.0_rfreal
217  w = 0.0_rfreal
218  p = pmixtinput%prepRealVal4
219  ELSE
220  d = pmixtinput%prepRealVal5
221  u = 0.0_rfreal
222  v = 0.0_rfreal
223  w = 0.0_rfreal
224  p = pmixtinput%prepRealVal6
225  END IF ! x
226 
227  mw = pgv(gv_mixt_mol,indmol*icg)
228  cp = pgv(gv_mixt_cp ,indcp *icg)
229 
230  gc = mixtperf_r_m(mw)
231  g = mixtperf_g_cpr(cp,gc)
232 
233  pcv(cv_mixt_dens,icg) = d
234  pcv(cv_mixt_xmom,icg) = d*u
235  pcv(cv_mixt_ymom,icg) = d*v
236  pcv(cv_mixt_zmom,icg) = d*w
237  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
238  END DO ! icg
239 
240 ! ----- With backchamber -----------------------------------------------------
241 
242  CASE ( "ds_7p25_v6", "ds_20p0_v6", "ds_50p0_v6", "ds_125p0_v6" )
243  DO icg = 1,pgrid%nCellsTot
244  x = pgrid%cofg(xcoord,icg)
245  y = abs(pgrid%cofg(ycoord,icg))
246 
247  IF ( (x < pmixtinput%prepRealVal1) .AND. &
248  (y < pmixtinput%prepRealVal2) ) THEN
249  d = pmixtinput%prepRealVal3
250  u = pmixtinput%prepRealVal4
251  v = 0.0_rfreal
252  w = 0.0_rfreal
253  p = pmixtinput%prepRealVal5
254  ELSE
255  d = pmixtinput%prepRealVal6
256  u = 0.0_rfreal
257  v = 0.0_rfreal
258  w = 0.0_rfreal
259  p = pmixtinput%prepRealVal7
260  END IF ! x
261 
262  mw = pgv(gv_mixt_mol,indmol*icg)
263  cp = pgv(gv_mixt_cp ,indcp *icg)
264 
265  gc = mixtperf_r_m(mw)
266  g = mixtperf_g_cpr(cp,gc)
267 
268  pcv(cv_mixt_dens,icg) = d
269  pcv(cv_mixt_xmom,icg) = d*u
270  pcv(cv_mixt_ymom,icg) = d*v
271  pcv(cv_mixt_zmom,icg) = d*w
272  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
273  END DO ! icg
274 
275 ! ------------------------------------------------------------------------------
276 ! hexahedral channel mesh
277 ! ------------------------------------------------------------------------------
278 
279  CASE ( "vortexNSCBC" )
280 
281 ! a = speed of sound
282 ! C = -0.0005_RFREAL*a
283 ! Rc = 0.15_RFREAL
284 ! u0 = 1.1_RFREAL*a
285 ! pinf = pMixtInput%prepRealVal3
286 
287  DO icg = 1,pgrid%nCellsTot
288  x = pgrid%cofg(xcoord,icg)
289  y = pgrid%cofg(ycoord,icg)
290 
291 ! radius = (x*x+y*y)**(0.5_RFREAL)
292 ! psi = C*exp(-(x*x+y*y)/(2.0_RFREAL*Rc*Rc))
293 
294 ! IF ( radius < Rc ) THEN
295 ! d = pMixtInput%prepRealVal2
296 ! u = u0 - x*psi/(d*Rc*Rc)
297 ! v = y*psi/(d*Rc*Rc)
298 ! w = 0.0_RFREAL
299 ! p = pinf + d*C*psi/(Rc*Rc)
300 ! ELSE
301 ! d = pMixtInput%prepRealVal2
302 ! u = u0
303 ! v = 0.0_RFREAL
304 ! w = 0.0_RFREAL
305 ! p = pinf
306 ! END IF ! radius
307 
308  mw = pgv(gv_mixt_mol,indmol*icg)
309  cp = pgv(gv_mixt_cp ,indcp *icg)
310 
311  gc = mixtperf_r_m(mw)
312  g = mixtperf_g_cpr(cp,gc)
313 
314  pcv(cv_mixt_dens,icg) = d
315  pcv(cv_mixt_xmom,icg) = d*u
316  pcv(cv_mixt_ymom,icg) = d*v
317  pcv(cv_mixt_zmom,icg) = d*w
318  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
319  END DO ! icg
320 
321 ! ------------------------------------------------------------------------------
322 ! Generic compressible gravity current
323 ! ------------------------------------------------------------------------------
324 
325  CASE ( "gcgc" )
326  DO icg = 1,pgrid%nCellsTot
327  x = pgrid%cofg(xcoord,icg)
328  y = pgrid%cofg(ycoord,icg)
329 
330  IF ( (x < pmixtinput%prepRealVal1) .AND. &
331  (y < pmixtinput%prepRealVal2) ) THEN
332  d = pmixtinput%prepRealVal3
333  u = 0.0_rfreal
334  v = 0.0_rfreal
335  w = 0.0_rfreal
336  p = pmixtinput%prepRealVal4
337  ELSE
338  d = pmixtinput%prepRealVal5
339  u = 0.0_rfreal
340  v = 0.0_rfreal
341  w = 0.0_rfreal
342  p = pmixtinput%prepRealVal4
343  END IF ! x
344 
345  mw = pgv(gv_mixt_mol,indmol*icg)
346  cp = pgv(gv_mixt_cp ,indcp *icg)
347 
348  gc = mixtperf_r_m(mw)
349  g = mixtperf_g_cpr(cp,gc)
350 
351  pcv(cv_mixt_dens,icg) = d
352  pcv(cv_mixt_xmom,icg) = d*u
353  pcv(cv_mixt_ymom,icg) = d*v
354  pcv(cv_mixt_zmom,icg) = d*w
355  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
356  END DO ! icg
357 
358 ! ------------------------------------------------------------------------------
359 ! Generic multiphase wall jet
360 ! ------------------------------------------------------------------------------
361 
362  CASE ( "gmpjet" )
363  DO icg = 1,pgrid%nCellsTot
364  x = pgrid%cofg(xcoord,icg)
365  y = pgrid%cofg(ycoord,icg)
366  z = pgrid%cofg(zcoord,icg)
367 
368  IF ( ( x < 0.00_rfreal) .AND. &
369  (abs(y) < 0.55_rfreal) .AND. &
370  (abs(z) < 0.55_rfreal) ) THEN
371  d = 34.7418238572_rfreal
372  u = 0.0_rfreal
373  v = 0.0_rfreal
374  w = 0.0_rfreal
375  p = 5.0e+6_rfreal
376  ELSE
377  d = 1.20903522664_rfreal
378  u = 0.0_rfreal
379  v = 0.0_rfreal
380  w = 0.0_rfreal
381  p = 1.0e+5_rfreal
382  END IF ! x
383 
384  mw = pgv(gv_mixt_mol,indmol*icg)
385  cp = pgv(gv_mixt_cp ,indcp *icg)
386 
387  gc = mixtperf_r_m(mw)
388  g = mixtperf_g_cpr(cp,gc)
389 
390  pcv(cv_mixt_dens,icg) = d
391  pcv(cv_mixt_xmom,icg) = d*u
392  pcv(cv_mixt_ymom,icg) = d*v
393  pcv(cv_mixt_zmom,icg) = d*w
394  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
395  END DO ! icg
396 
397 ! ------------------------------------------------------------------------------
398 ! Gradient test
399 ! ------------------------------------------------------------------------------
400 
401 ! ----- Linear function --------------------------------------------------------
402 
403  CASE ( "gtlin" )
404  xmin = minval(pgrid%cofg(xcoord,1:pgrid%nCellsTot))
405  ymin = minval(pgrid%cofg(ycoord,1:pgrid%nCellsTot))
406  zmin = minval(pgrid%cofg(zcoord,1:pgrid%nCellsTot))
407 
408  CALL rflu_setexactflowlinear(xmin,ymin,zmin,1,dmin,gx,gy,gz)
409  CALL rflu_setexactflowlinear(xmin,ymin,zmin,5,pmin,gx,gy,gz)
410 
411  IF ( dmin < 0.0_rfreal ) THEN
412  doffs = -2.0_rfreal*dmin
413  ELSE
414  doffs = 0.0_rfreal
415  END IF ! dMin
416 
417  IF ( pmin < 0.0_rfreal ) THEN
418  poffs = -2.0_rfreal*pmin
419  ELSE
420  poffs = 0.0_rfreal
421  END IF ! pMin
422 
423  DO icg = 1,pgrid%nCellsTot
424  x = pgrid%cofg(xcoord,icg)
425  y = pgrid%cofg(ycoord,icg)
426  z = pgrid%cofg(zcoord,icg)
427 
428  mw = pgv(gv_mixt_mol,indmol*icg)
429  cp = pgv(gv_mixt_cp ,indcp *icg)
430 
431  gc = mixtperf_r_m(mw)
432  g = mixtperf_g_cpr(cp,gc)
433 
434  CALL rflu_setexactflowlinear(x,y,z,1,d,gx,gy,gz)
435  CALL rflu_setexactflowlinear(x,y,z,2,u,gx,gy,gz)
436  CALL rflu_setexactflowlinear(x,y,z,3,v,gx,gy,gz)
437  CALL rflu_setexactflowlinear(x,y,z,4,w,gx,gy,gz)
438  CALL rflu_setexactflowlinear(x,y,z,5,p,gx,gy,gz)
439 
440  d = d + doffs
441  p = p + poffs
442 
443  pcv(cv_mixt_dens,icg) = d
444  pcv(cv_mixt_xmom,icg) = d*u
445  pcv(cv_mixt_ymom,icg) = d*v
446  pcv(cv_mixt_zmom,icg) = d*w
447  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
448  END DO ! icg
449 
450 ! ----- Trigonometric function -------------------------------------------------
451 
452  CASE ( "gttri" )
453  nx = pregion%mixtInput%prepRealVal1
454  ny = pregion%mixtInput%prepRealVal2
455  nz = pregion%mixtInput%prepRealVal3
456 
457  doffs = 2.0_rfreal ! Assume amplitude is unity
458  poffs = 2.0_rfreal ! Assume amplitude is unity
459 
460  DO icg = 1,pgrid%nCellsTot
461  x = pgrid%cofg(xcoord,icg)
462  y = pgrid%cofg(ycoord,icg)
463  z = pgrid%cofg(zcoord,icg)
464 
465  mw = pgv(gv_mixt_mol,indmol*icg)
466  cp = pgv(gv_mixt_cp ,indcp *icg)
467 
468  gc = mixtperf_r_m(mw)
469  g = mixtperf_g_cpr(cp,gc)
470 
471  CALL rflu_setexactflowtrig(global,nx,ny,nz,x,y,z,1,d,gx,gy,gz)
472  CALL rflu_setexactflowtrig(global,nx,ny,nz,x,y,z,2,u,gx,gy,gz)
473  CALL rflu_setexactflowtrig(global,nx,ny,nz,x,y,z,3,v,gx,gy,gz)
474  CALL rflu_setexactflowtrig(global,nx,ny,nz,x,y,z,4,w,gx,gy,gz)
475  CALL rflu_setexactflowtrig(global,nx,ny,nz,x,y,z,5,p,gx,gy,gz)
476 
477  d = d + doffs
478  p = p + poffs
479 
480  pcv(cv_mixt_dens,icg) = d
481  pcv(cv_mixt_xmom,icg) = d*u
482  pcv(cv_mixt_ymom,icg) = d*v
483  pcv(cv_mixt_zmom,icg) = d*w
484  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
485  END DO ! icg
486 
487 ! ------------------------------------------------------------------------------
488 ! Jet Flow
489 ! ------------------------------------------------------------------------------
490 
491  CASE ( "jet" )
492  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
493  CALL errorstop(global,err_gasmodel_invalid,__line__, &
494  'Case initialization only valid with gas-liq model.')
495  END IF ! pRegion%mixtInput%gasModel
496 
497  IF ( pregion%specInput%nSpecies /= 2 ) THEN
498  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
499  pregion%specInput%nSpecies
500  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
501  trim(errorstring))
502  END IF ! pRegion%specInput%nSpecies
503 
504  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
505  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
506  rgas)
507  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
508  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
509  rvap)
510  cvl = global%refCvLiq
511 
512  DO icg = 1,pgrid%nCellsTot
513  x = pgrid%cofg(xcoord,icg)
514  y = pgrid%cofg(ycoord,icg)
515 
516  yp = 8.9e-05_rfreal + 0.0175_rfreal*(x - 2.5e-04_rfreal)
517 
518  IF ( (y <= 8.9e-05_rfreal) .AND. (y >= -8.9e-05_rfreal) .AND. &
519  (x >= 0.0e+00_rfreal) .AND. (x <= 2.5e-04_rfreal) ) THEN
520 
521  rl = 880.0_rfreal
522  pl = 1.0_rfreal
523  rg = 11.69_rfreal
524  pg = 0.0_rfreal
525  rv = 0.722_rfreal
526  pv = 0.0_rfreal
527 
528  d = rl*pl + rg*pg + rv*pv
529  h = 8.9e-05_rfreal
530  u = 300.0*(1.0_rfreal - (y*y)/(h*h))
531 
532  IF ( abs(y) < 5.5e-05_rfreal ) THEN
533  u = 300.0_rfreal
534  END IF ! ABS(y)
535 
536  v = 0.0_rfreal
537  w = 0.0_rfreal
538  t = 288.0_rfreal
539  ELSE
540  rl = 880.0_rfreal
541  pl = 0.0_rfreal
542  rg = 11.69_rfreal
543  pg = 1.0_rfreal
544  rv = 0.722_rfreal
545  pv = 0.0_rfreal
546 
547  d = rl*pl + rg*pg + rv*pv
548  u = 0.0_rfreal
549  v = 0.0_rfreal
550  w = 0.0_rfreal
551  t = 288.0_rfreal
552  END IF ! y
553 
554  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
555  vm2 = u*u+v*v+w*w
556 
557  pcv(cv_mixt_dens,icg) = d
558  pcv(cv_mixt_xmom,icg) = d*u
559  pcv(cv_mixt_ymom,icg) = d*v
560  pcv(cv_mixt_zmom,icg) = d*w
561  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
562  END DO ! icg
563 
564 ! ------------------------------------------------------------------------------
565 ! Kieffer jet. NOTE configuration 1 has jet axis along z, configuration 2
566 ! has jet axis along x.
567 ! ------------------------------------------------------------------------------
568 
569  CASE ( "kjet" )
570  DO icg = 1,pgrid%nCellsTot
571  z = pgrid%cofg(zcoord,icg)
572 
573  IF ( z < -0.0254_rfreal ) THEN
574  d = 7.96_rfreal
575  u = 0.0_rfreal
576  v = 0.0_rfreal
577  w = 0.0_rfreal
578  p = 7.25e+5_rfreal
579  ELSE
580  d = 1.16871466381_rfreal
581  u = 0.0_rfreal
582  v = 0.0_rfreal
583  w = 0.0_rfreal
584  p = 1.0e+5_rfreal
585  END IF ! z
586 
587  mw = pgv(gv_mixt_mol,indmol*icg)
588  cp = pgv(gv_mixt_cp ,indcp *icg)
589 
590  gc = mixtperf_r_m(mw)
591  g = mixtperf_g_cpr(cp,gc)
592 
593  pcv(cv_mixt_dens,icg) = d
594  pcv(cv_mixt_xmom,icg) = d*u
595  pcv(cv_mixt_ymom,icg) = d*v
596  pcv(cv_mixt_zmom,icg) = d*w
597  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
598  END DO ! icg
599 
600  CASE ( "kjet2","kjet2v3","kjet2v4","kjet2v5" )
601  DO icg = 1,pgrid%nCellsTot
602  x = pgrid%cofg(xcoord,icg)
603 
604  IF ( x < -0.0127_rfreal ) THEN
605  d = 7.96_rfreal
606  u = 0.0_rfreal
607  v = 0.0_rfreal
608  w = 0.0_rfreal
609  p = 7.25e+5_rfreal
610  ELSE
611  d = 1.1220002356_rfreal
612  u = 0.0_rfreal
613  v = 0.0_rfreal
614  w = 0.0_rfreal
615  p = 1.0e+5_rfreal
616  END IF ! x
617 
618  mw = pgv(gv_mixt_mol,indmol*icg)
619  cp = pgv(gv_mixt_cp ,indcp *icg)
620 
621  gc = mixtperf_r_m(mw)
622  g = mixtperf_g_cpr(cp,gc)
623 
624  pcv(cv_mixt_dens,icg) = d
625  pcv(cv_mixt_xmom,icg) = d*u
626  pcv(cv_mixt_ymom,icg) = d*v
627  pcv(cv_mixt_zmom,icg) = d*w
628  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
629  END DO ! icg
630 
631  CASE ( "kjet2v3mp","kjet2v4mp","kjet2v5mp" )
632  DO icg = 1,pgrid%nCellsTot
633  x = pgrid%cofg(xcoord,icg)
634 
635  IF ( x < pmixtinput%prepRealVal1 ) THEN
636  d = pmixtinput%prepRealVal2
637  u = pmixtinput%prepRealVal3
638  v = 0.0_rfreal
639  w = 0.0_rfreal
640  p = pmixtinput%prepRealVal4
641  ELSE
642  d = pmixtinput%prepRealVal5
643  u = pmixtinput%prepRealVal6
644  v = 0.0_rfreal
645  w = 0.0_rfreal
646  p = pmixtinput%prepRealVal7
647  END IF ! x
648 
649  mw = pgv(gv_mixt_mol,indmol*icg)
650  cp = pgv(gv_mixt_cp ,indcp *icg)
651 
652  gc = mixtperf_r_m(mw)
653  g = mixtperf_g_cpr(cp,gc)
654 
655  pcv(cv_mixt_dens,icg) = d
656  pcv(cv_mixt_xmom,icg) = d*u
657  pcv(cv_mixt_ymom,icg) = d*v
658  pcv(cv_mixt_zmom,icg) = d*w
659  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
660  END DO ! icg
661 
662 ! ------------------------------------------------------------------------------
663 ! Multiphase Shocktube: Water-Air
664 ! ------------------------------------------------------------------------------
665 
666  CASE ( "MShock_H2O_Air001" )
667  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
668  CALL errorstop(global,err_gasmodel_invalid,__line__, &
669  'Case initialization only valid with gas-liq model.')
670  END IF ! pRegion%mixtInput%gasModel
671 
672  IF ( pregion%specInput%nSpecies /= 2 ) THEN
673  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
674  pregion%specInput%nSpecies
675  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
676  trim(errorstring))
677  END IF ! pRegion%specInput%nSpecies
678 
679  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
680  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
681  rgas)
682  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
683  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
684  rvap)
685  cvl = global%refCvLiq
686 
687  DO icg = 1,pgrid%nCellsTot
688  x = pgrid%cofg(xcoord,icg)
689 
690  IF ( x < 0.7_rfreal ) THEN
691  rl = 1500.0_rfreal
692  pl = 1.0_rfreal
693  rg = 0.0_rfreal
694  pg = 0.0_rfreal
695  rv = 0.0_rfreal
696  pv = 0.0_rfreal
697 
698  d = rl*pl + rg*pg + rv*pv
699  u = 0.0_rfreal
700  v = 0.0_rfreal
701  w = 0.0_rfreal
702  t = 515.0_rfreal
703  ELSE
704  rl = 0.0_rfreal
705  pl = 0.0_rfreal
706  rg = 150.0_rfreal
707  pg = 1.0_rfreal
708  rv = 0.0_rfreal
709  pv = 0.0_rfreal
710 
711  d = rl*pl + rg*pg + rv*pv
712  u = 0.0_rfreal
713  v = 0.0_rfreal
714  w = 0.0_rfreal
715  t = 2.25_rfreal
716  END IF ! x
717 
718  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
719  vm2 = u*u+v*v+w*w
720 
721  pcv(cv_mixt_dens,icg) = d
722  pcv(cv_mixt_xmom,icg) = d*u
723  pcv(cv_mixt_ymom,icg) = d*v
724  pcv(cv_mixt_zmom,icg) = d*w
725  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
726  END DO ! icg
727 
728  CASE ( "MShock_H2O_Air002" )
729  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
730  CALL errorstop(global,err_gasmodel_invalid,__line__, &
731  'Case initialization only valid with gas-liq model.')
732  END IF ! pRegion%mixtInput%gasModel
733 
734  IF ( pregion%specInput%nSpecies /= 2 ) THEN
735  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
736  pregion%specInput%nSpecies
737  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
738  trim(errorstring))
739  END IF ! pRegion%specInput%nSpecies
740 
741  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
742  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
743  rgas)
744  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
745  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
746  rvap)
747  cvl = global%refCvLiq
748 
749  DO icg = 1,pgrid%nCellsTot
750  x = pgrid%cofg(xcoord,icg)
751 
752  IF ( x < 0.7_rfreal ) THEN
753  rl = 1500.0_rfreal
754  pl = 1.0_rfreal
755  rg = 0.0_rfreal
756  pg = 0.0_rfreal
757  rv = 0.0_rfreal
758  pv = 0.0_rfreal
759 
760  d = rl*pl + rg*pg + rv*pv
761  u = 0.0_rfreal
762  v = 0.0_rfreal
763  w = 0.0_rfreal
764  t = 515.0_rfreal
765  ELSE
766  rl = 0.0_rfreal
767  pl = 0.0_rfreal
768  rg = 50.0_rfreal
769  pg = 1.0_rfreal
770  rv = 0.0_rfreal
771  pv = 0.0_rfreal
772 
773  d = rl*pl + rg*pg + rv*pv
774  u = 0.0_rfreal
775  v = 0.0_rfreal
776  w = 0.0_rfreal
777  t = 6.73_rfreal
778  END IF ! x
779 
780  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
781  vm2 = u*u+v*v+w*w
782 
783  pcv(cv_mixt_dens,icg) = d
784  pcv(cv_mixt_xmom,icg) = d*u
785  pcv(cv_mixt_ymom,icg) = d*v
786  pcv(cv_mixt_zmom,icg) = d*w
787  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
788  END DO ! icg
789 
790 ! ------------------------------------------------------------------------------
791 ! Multiphase Shocktube: Air-Air-He
792 ! ------------------------------------------------------------------------------
793 
794  CASE ( "MShock_Air_Air_He001" )
795  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
796  CALL errorstop(global,err_gasmodel_invalid,__line__, &
797  'Case initialization only valid with gas-liq model.')
798  END IF ! pRegion%mixtInput%gasModel
799 
800  IF ( pregion%specInput%nSpecies /= 2 ) THEN
801  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
802  pregion%specInput%nSpecies
803  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
804  trim(errorstring))
805  END IF ! pRegion%specInput%nSpecies
806 
807  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
808  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
809  rgas)
810  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
811  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
812  rvap)
813  cvl = global%refCvLiq
814 
815  DO icg = 1,pgrid%nCellsTot
816  x = pgrid%cofg(xcoord,icg)
817 
818  IF ( x < 0.25_rfreal ) THEN
819  rl = 0.0_rfreal
820  pl = 0.0_rfreal
821  rg = 1.7017_rfreal
822  pg = 1.0_rfreal
823  rv = 0.0_rfreal
824  pv = 0.0_rfreal
825 
826  d = rl*pl + rg*pg + rv*pv
827  u = 98.956_rfreal
828  v = 0.0_rfreal
829  w = 0.0_rfreal
830  t = 307.13_rfreal
831  ELSE IF ( (x > 0.25_rfreal) .AND. (x < 0.5_rfreal) ) THEN
832  rl = 0.0_rfreal
833  pl = 0.0_rfreal
834  rg = 1.2763
835  pg = 1.0_rfreal
836  rv = 0.0_rfreal
837  pv = 0.0_rfreal
838 
839  d = rl*pl + rg*pg + rv*pv
840  u = 0.0_rfreal
841  v = 0.0_rfreal
842  w = 0.0_rfreal
843  t = 273.0_rfreal
844  ELSE
845  rl = 0.0_rfreal
846  pl = 0.0_rfreal
847  rg = 0.0_rfreal
848  pg = 0.0_rfreal
849  rv = 0.1760_rfreal
850  pv = 1.0_rfreal
851 
852  d = rl*pl + rg*pg + rv*pv
853  u = 0.0_rfreal
854  v = 0.0_rfreal
855  w = 0.0_rfreal
856  t = 273.5588_rfreal
857  END IF ! x
858 
859  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
860  vm2 = u*u+v*v+w*w
861 
862  pcv(cv_mixt_dens,icg) = d
863  pcv(cv_mixt_xmom,icg) = d*u
864  pcv(cv_mixt_ymom,icg) = d*v
865  pcv(cv_mixt_zmom,icg) = d*w
866  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
867  END DO ! icg
868 
869 ! ------------------------------------------------------------------------------
870 ! Moving shock
871 ! ------------------------------------------------------------------------------
872 
873  CASE ( "mvsh" )
874  DO icg = 1,pgrid%nCellsTot
875  x = pgrid%cofg(xcoord,icg)
876 
877  IF ( x < pmixtinput%prepRealVal1 ) THEN
878  d = pmixtinput%prepRealVal2
879  u = pmixtinput%prepRealVal3
880  v = 0.0_rfreal
881  w = 0.0_rfreal
882  p = pmixtinput%prepRealVal4
883  ELSE
884  d = pmixtinput%prepRealVal5
885  u = 0.0_rfreal
886  v = 0.0_rfreal
887  w = 0.0_rfreal
888  p = pmixtinput%prepRealVal6
889  END IF ! z
890 
891  mw = pgv(gv_mixt_mol,indmol*icg)
892  cp = pgv(gv_mixt_cp ,indcp *icg)
893 
894  gc = mixtperf_r_m(mw)
895  g = mixtperf_g_cpr(cp,gc)
896 
897  pcv(cv_mixt_dens,icg) = d
898  pcv(cv_mixt_xmom,icg) = d*u
899  pcv(cv_mixt_ymom,icg) = d*v
900  pcv(cv_mixt_zmom,icg) = d*w
901  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
902  END DO ! icg
903 
904 ! ------------------------------------------------------------------------------
905 ! NSCBC farfield test
906 ! ------------------------------------------------------------------------------
907 
908  CASE ( "farf" )
909  l_l = 2.0_rfreal
910  a_c = 0.0001_rfreal
911  a_c = 0.01_rfreal
912  uo_c = 0.1_rfreal
913  uo_c = 0.0_rfreal
914 
915  r1 = 0.006125_rfreal
916  r2 = 0.091875_rfreal
917 
918  l = l_l*pmixtinput%prepRealVal1
919  c = pmixtinput%prepRealVal2
920  a = a_c*c
921  uo = uo_c*c
922 
923  ro = pmixtinput%prepRealVal3
924  po = pmixtinput%prepRealVal4
925 
926  DO icg = 1,pgrid%nCellsTot
927  x = pgrid%cofg(xcoord,icg)
928  y = pgrid%cofg(ycoord,icg)
929 
930  radius = sqrt(x*x+y*y)
931 
932  mw = pgv(gv_mixt_mol,indmol*icg)
933  cp = pgv(gv_mixt_cp ,indcp *icg)
934 
935  gc = mixtperf_r_m(mw)
936  g = mixtperf_g_cpr(cp,gc)
937 
938  um= uo + a*(sin(global%pi*(radius-r1)/(r2-r1)))**9.0_rfreal
939  u = um*cos(atan2(y,x))
940  p = po + (um - uo)*(ro*c)
941  d = ro*(p/po)**(1.0_rfreal/g)
942  v = um*sin(atan2(y,x))
943  w = 0.0_rfreal
944 
945  pcv(cv_mixt_dens,icg) = d
946  pcv(cv_mixt_xmom,icg) = d*u
947  pcv(cv_mixt_ymom,icg) = d*v
948  pcv(cv_mixt_zmom,icg) = d*w
949  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
950  END DO ! icg
951 
952 ! ------------------------------------------------------------------------------
953 ! NSCBC tests
954 ! ------------------------------------------------------------------------------
955 
956  CASE ( "nscbc1" )
957  l_l = 2.0_rfreal
958  a_c = 0.0001_rfreal
959  uo_c = 0.1_rfreal
960 
961  l = l_l*pmixtinput%prepRealVal1
962  c = pmixtinput%prepRealVal2
963  a = a_c*c
964  uo = uo_c*c
965 
966  ro = pmixtinput%prepRealVal3
967  po = pmixtinput%prepRealVal4
968 
969  DO icg = 1,pgrid%nCellsTot
970  x = pgrid%cofg(xcoord,icg)
971 
972  mw = pgv(gv_mixt_mol,indmol*icg)
973  cp = pgv(gv_mixt_cp ,indcp *icg)
974 
975  gc = mixtperf_r_m(mw)
976  g = mixtperf_g_cpr(cp,gc)
977 
978  u = uo + a*(sin(global%pi*x/l))**9.0_rfreal
979  p = po + (u - uo)*(ro*c)
980  d = ro*(p/po)**(1.0_rfreal/g)
981  v = 0.0_rfreal
982  w = 0.0_rfreal
983 
984  pcv(cv_mixt_dens,icg) = d
985  pcv(cv_mixt_xmom,icg) = d*u
986  pcv(cv_mixt_ymom,icg) = d*v
987  pcv(cv_mixt_zmom,icg) = d*w
988  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
989  END DO ! icg
990 
991  CASE ( "nscbc2" )
992  l_l = 2.0_rfreal
993  a_c = 0.01_rfreal
994  uo_c = 0.1_rfreal
995 
996  l = l_l*pmixtinput%prepRealVal1
997  c = pmixtinput%prepRealVal2
998  a = a_c*c
999  uo = uo_c*c
1000 
1001  ro = pmixtinput%prepRealVal3
1002  po = pmixtinput%prepRealVal4
1003 
1004  DO icg = 1,pgrid%nCellsTot
1005  x = pgrid%cofg(xcoord,icg)
1006 
1007  mw = pgv(gv_mixt_mol,indmol*icg)
1008  cp = pgv(gv_mixt_cp ,indcp *icg)
1009 
1010  gc = mixtperf_r_m(mw)
1011  g = mixtperf_g_cpr(cp,gc)
1012 
1013  u = uo
1014  p = po + ro*c*a*((sin(global%pi*x/l))**9.0_rfreal)
1015  d = ro*(p/po)**(1.0_rfreal/g)
1016  v = 0.0_rfreal
1017  w = 0.0_rfreal
1018 
1019  pcv(cv_mixt_dens,icg) = d
1020  pcv(cv_mixt_xmom,icg) = d*u
1021  pcv(cv_mixt_ymom,icg) = d*v
1022  pcv(cv_mixt_zmom,icg) = d*w
1023  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1024  END DO ! icg
1025 
1026  CASE ( "nscbc3" )
1027  l_l = 2.0_rfreal
1028  a_c = 0.0001_rfreal
1029  uo_c = 0.1_rfreal
1030 
1031  l = l_l*pmixtinput%prepRealVal1
1032  c = pmixtinput%prepRealVal2
1033  a = a_c*c
1034  uo = uo_c*c
1035 
1036  ro = pmixtinput%prepRealVal3
1037  po = pmixtinput%prepRealVal4
1038 
1039  DO icg = 1,pgrid%nCellsTot
1040  y = pgrid%cofg(ycoord,icg)
1041 
1042  mw = pgv(gv_mixt_mol,indmol*icg)
1043  cp = pgv(gv_mixt_cp ,indcp *icg)
1044 
1045  gc = mixtperf_r_m(mw)
1046  g = mixtperf_g_cpr(cp,gc)
1047 
1048  v = uo + a*(sin(global%pi*y/l))**9.0_rfreal
1049  p = po + (v - uo)*(ro*c)
1050  d = ro*(p/po)**(1.0_rfreal/g)
1051  u = 0.0_rfreal
1052  w = 0.0_rfreal
1053 
1054  pcv(cv_mixt_dens,icg) = d
1055  pcv(cv_mixt_xmom,icg) = d*u
1056  pcv(cv_mixt_ymom,icg) = d*v
1057  pcv(cv_mixt_zmom,icg) = d*w
1058  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1059  END DO ! icg
1060 
1061  CASE ( "nscbc4" )
1062  l_l = 2.0_rfreal
1063  a_c = 0.01_rfreal
1064  uo_c = 0.1_rfreal
1065 
1066  l = l_l*pmixtinput%prepRealVal1
1067  c = pmixtinput%prepRealVal2
1068  a = a_c*c
1069  uo = uo_c*c
1070 
1071  ro = pmixtinput%prepRealVal3
1072  po = pmixtinput%prepRealVal4
1073 
1074  DO icg = 1,pgrid%nCellsTot
1075  y = pgrid%cofg(ycoord,icg)
1076 
1077  mw = pgv(gv_mixt_mol,indmol*icg)
1078  cp = pgv(gv_mixt_cp ,indcp *icg)
1079 
1080  gc = mixtperf_r_m(mw)
1081  g = mixtperf_g_cpr(cp,gc)
1082 
1083  u = 0.0_rfreal
1084  p = po + ro*c*a*((sin(global%pi*y/l))**9.0_rfreal)
1085  d = ro*(p/po)**(1.0_rfreal/g)
1086  v = uo
1087  w = 0.0_rfreal
1088 
1089  pcv(cv_mixt_dens,icg) = d
1090  pcv(cv_mixt_xmom,icg) = d*u
1091  pcv(cv_mixt_ymom,icg) = d*v
1092  pcv(cv_mixt_zmom,icg) = d*w
1093  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1094  END DO ! icg
1095 
1096  CASE ( "nscbc5" )
1097  l_l = 2.0_rfreal
1098  a_c = 0.0001_rfreal
1099  uo_c = 0.1_rfreal
1100 
1101  l = l_l*pmixtinput%prepRealVal1
1102  c = pmixtinput%prepRealVal2
1103  a = a_c*c
1104  uo = uo_c*c
1105 
1106  ro = pmixtinput%prepRealVal3
1107  po = pmixtinput%prepRealVal4
1108 
1109  DO icg = 1,pgrid%nCellsTot
1110  x = pgrid%cofg(xcoord,icg)
1111  y = pgrid%cofg(ycoord,icg)
1112 
1113  xx = ((x*x+y*y)**0.5)*cos(atan2(y,x)-global%pi/4.0_rfreal)
1114 
1115  mw = pgv(gv_mixt_mol,indmol*icg)
1116  cp = pgv(gv_mixt_cp ,indcp *icg)
1117 
1118  gc = mixtperf_r_m(mw)
1119  g = mixtperf_g_cpr(cp,gc)
1120 
1121  um = uo + a*(sin(global%pi*xx/l))**9.0_rfreal
1122 
1123  u = um*cos(global%pi/4.0_rfreal)
1124  p = po + (um - uo)*(ro*c)
1125  d = ro*(p/po)**(1.0_rfreal/g)
1126  v = um*sin(global%pi/4.0_rfreal)
1127  w = 0.0_rfreal
1128 
1129  pcv(cv_mixt_dens,icg) = d
1130  pcv(cv_mixt_xmom,icg) = d*u
1131  pcv(cv_mixt_ymom,icg) = d*v
1132  pcv(cv_mixt_zmom,icg) = d*w
1133  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1134  END DO ! icg
1135 
1136  CASE ( "nscbc6" )
1137  l_l = 2.0_rfreal
1138  a_c = 0.01_rfreal
1139  uo_c = 0.1_rfreal
1140 
1141  l = l_l*pmixtinput%prepRealVal1
1142  c = pmixtinput%prepRealVal2
1143  a = a_c*c
1144  uo = uo_c*c
1145 
1146  ro = pmixtinput%prepRealVal3
1147  po = pmixtinput%prepRealVal4
1148 
1149  DO icg = 1,pgrid%nCellsTot
1150  x = pgrid%cofg(xcoord,icg)
1151  y = pgrid%cofg(ycoord,icg)
1152 
1153  xx = ((x*x+y*y)**0.5)*cos(atan2(y,x)-global%pi/4.0_rfreal)
1154 
1155  mw = pgv(gv_mixt_mol,indmol*icg)
1156  cp = pgv(gv_mixt_cp ,indcp *icg)
1157 
1158  gc = mixtperf_r_m(mw)
1159  g = mixtperf_g_cpr(cp,gc)
1160 
1161  um = uo + a*(sin(global%pi*xx/l))**9.0_rfreal
1162 
1163  u = uo*cos(global%pi/4.0_rfreal)
1164  p = po + ro*c*a*((sin(global%pi*xx/l))**9.0_rfreal)
1165  d = ro*(p/po)**(1.0_rfreal/g)
1166  v = uo*sin(global%pi/4.0_rfreal)
1167  w = 0.0_rfreal
1168 
1169  pcv(cv_mixt_dens,icg) = d
1170  pcv(cv_mixt_xmom,icg) = d*u
1171  pcv(cv_mixt_ymom,icg) = d*v
1172  pcv(cv_mixt_zmom,icg) = d*w
1173  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1174  END DO ! icg
1175 
1176  CASE ( "nscbc7" )
1177  l_l = 2.0_rfreal
1178  a_c = 0.01_rfreal
1179  uo_c = 0.1_rfreal
1180 
1181  l = l_l*pmixtinput%prepRealVal1
1182  c = pmixtinput%prepRealVal2
1183  a = a_c*c
1184  uo = uo_c*c
1185 
1186  ro = pmixtinput%prepRealVal3
1187  po = pmixtinput%prepRealVal4
1188  xo = pmixtinput%prepRealVal5
1189 
1190  DO icg = 1,pgrid%nCellsTot
1191  x = pgrid%cofg(xcoord,icg)
1192 
1193  mw = pgv(gv_mixt_mol,indmol*icg)
1194  cp = pgv(gv_mixt_cp ,indcp *icg)
1195 
1196  gc = mixtperf_r_m(mw)
1197  g = mixtperf_g_cpr(cp,gc)
1198 
1199  x = x + xo
1200 
1201  IF ( x < l ) THEN
1202  u = uo
1203  p = po + ro*c*a*((sin(global%pi*x/l))**9.0_rfreal)
1204  d = ro*(p/po)**(1.0_rfreal/g)
1205  v = 0.0_rfreal
1206  w = 0.0_rfreal
1207  ELSE
1208  u = uo
1209  p = po
1210  d = ro
1211  v = 0.0_rfreal
1212  w = 0.0_rfreal
1213  END IF ! x
1214 
1215  pcv(cv_mixt_dens,icg) = d
1216  pcv(cv_mixt_xmom,icg) = d*u
1217  pcv(cv_mixt_ymom,icg) = d*v
1218  pcv(cv_mixt_zmom,icg) = d*w
1219  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1220  END DO ! icg
1221 
1222 ! ------------------------------------------------------------------------------
1223 ! Shock-tube testing on nscbc
1224 ! ------------------------------------------------------------------------------
1225 
1226  CASE ( "nscbc8" )
1227  DO icg = 1,pgrid%nCellsTot
1228  x = pgrid%cofg(xcoord,icg)
1229 
1230  IF ( x < pmixtinput%prepRealVal1 ) THEN
1231  d = pmixtinput%prepRealVal2
1232  u = 0.0_rfreal
1233  v = 0.0_rfreal
1234  w = 0.0_rfreal
1235  p = pmixtinput%prepRealVal3
1236  ELSE
1237  d = pmixtinput%prepRealVal4
1238  u = 0.0_rfreal
1239  v = 0.0_rfreal
1240  w = 0.0_rfreal
1241  p = pmixtinput%prepRealVal5
1242  END IF ! x
1243 
1244  mw = pgv(gv_mixt_mol,indmol*icg)
1245  cp = pgv(gv_mixt_cp ,indcp *icg)
1246 
1247  gc = mixtperf_r_m(mw)
1248  g = mixtperf_g_cpr(cp,gc)
1249 
1250  pcv(cv_mixt_dens,icg) = d
1251  pcv(cv_mixt_xmom,icg) = d*u
1252  pcv(cv_mixt_ymom,icg) = d*v
1253  pcv(cv_mixt_zmom,icg) = d*w
1254  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1255  END DO ! icg
1256 
1257 ! ------------------------------------------------------------------------------
1258 ! Nozzle cavity
1259 ! ------------------------------------------------------------------------------
1260 
1261  CASE ( "ncavity" )
1262  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
1263  CALL errorstop(global,err_gasmodel_invalid,__line__, &
1264  'Case initialization only valid with gas-liq model.')
1265  END IF ! pRegion%mixtInput%gasModel
1266 
1267  IF ( pregion%specInput%nSpecies /= 2 ) THEN
1268  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
1269  pregion%specInput%nSpecies
1270  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
1271  trim(errorstring))
1272  END IF ! pRegion%specInput%nSpecies
1273 
1274  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
1275  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
1276  rgas)
1277  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
1278  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
1279  rvap)
1280  cvl = global%refCvLiq
1281 
1282  DO icg = 1,pgrid%nCellsTot
1283  x = pgrid%cofg(xcoord,icg)
1284 
1285  IF ( x <= 3.0e-03_rfreal ) THEN
1286  rl = 882.655_rfreal
1287  pl = 1.0_rfreal
1288  rg = 24.14836_rfreal
1289  pg = 0.0_rfreal
1290  rv = 15.547_rfreal
1291  pv = 0.0_rfreal
1292 
1293  d = rl*pl + rg*pg + rv*pv
1294  u = 0.0_rfreal
1295  v = 0.0_rfreal
1296  w = 0.0_rfreal
1297  t = 293.0_rfreal
1298  ELSE
1299  rl = 880.0_rfreal
1300  pl = 0.0_rfreal
1301  rg = 24.14836_rfreal
1302  pg = 1.0_rfreal
1303  rv = 15.547_rfreal
1304  pv = 0.0_rfreal
1305 
1306  d = rl*pl + rg*pg + rv*pv
1307  u = 0.0_rfreal
1308  v = 0.0_rfreal
1309  w = 0.0_rfreal
1310  t = 293_rfreal
1311  END IF ! x
1312 
1313  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
1314  vm2 = u*u+v*v+w*w
1315 
1316  pcv(cv_mixt_dens,icg) = d
1317  pcv(cv_mixt_xmom,icg) = d*u
1318  pcv(cv_mixt_ymom,icg) = d*v
1319  pcv(cv_mixt_zmom,icg) = d*w
1320  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
1321  END DO ! icg
1322 
1323 ! ------------------------------------------------------------------------------
1324 ! Proudman-Culick flow. NOTE this problem is two-dimensional and assumed
1325 ! to lie in the x-y plane, and that the injection boundary is located at
1326 ! y = -height.
1327 ! ------------------------------------------------------------------------------
1328 
1329  CASE ( "onera_c0", "onera_c0_2d_100x50" )
1330  CALL rflu_getparamshardcodeproudman(dinc,minj,vinj,ptot)
1331 
1332  height = minval(pgrid%xyz(ycoord,1:pgrid%nVertTot))
1333 
1334  DO icg = 1,pgrid%nCellsTot
1335  x = pgrid%cofg(xcoord,icg)
1336  y = pgrid%cofg(ycoord,icg)
1337 
1338  CALL rflu_computeexactflowproudman(global,x,y,height,dinc,vinj, &
1339  ptot,d,u,v,w,p)
1340 
1341  mw = pgv(gv_mixt_mol,indmol*icg)
1342  cp = pgv(gv_mixt_cp ,indcp *icg)
1343 
1344  gc = mixtperf_r_m(mw)
1345  g = mixtperf_g_cpr(cp,gc)
1346 
1347  pcv(cv_mixt_dens,icg) = d
1348  pcv(cv_mixt_xmom,icg) = d*u
1349  pcv(cv_mixt_ymom,icg) = d*v
1350  pcv(cv_mixt_zmom,icg) = d*w
1351  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1352  END DO ! icg
1353 
1354 ! ------------------------------------------------------------------------------
1355 ! Culick flow
1356 ! ------------------------------------------------------------------------------
1357 
1358  CASE ( "onera_c0_3d" )
1359  DO icg = 1,pgrid%nCellsTot
1360  x = pgrid%cofg(xcoord,icg)
1361  y = pgrid%cofg(ycoord,icg)
1362  z = pgrid%cofg(zcoord,icg)
1363 
1364  CALL rflu_computeexactflowculick(global,x,y,z,d,u,v,w,p)
1365 
1366  mw = pgv(gv_mixt_mol,indmol*icg)
1367  cp = pgv(gv_mixt_cp ,indcp *icg)
1368 
1369  gc = mixtperf_r_m(mw)
1370  g = mixtperf_g_cpr(cp,gc)
1371 
1372  pcv(cv_mixt_dens,icg) = d
1373  pcv(cv_mixt_xmom,icg) = d*u
1374  pcv(cv_mixt_ymom,icg) = d*v
1375  pcv(cv_mixt_zmom,icg) = d*w
1376  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1377  END DO ! icg
1378 
1379 ! ------------------------------------------------------------------------------
1380 ! Pipe acoustics. NOTE the pipe is assumed to have the x-coordinate
1381 ! running down the axis.
1382 ! ------------------------------------------------------------------------------
1383 
1384  CASE ( "pipeacoust" )
1385  CALL rflu_getparamshardcodepacoust(ptot,atot)
1386 
1387  dtot = mixtperf_d_cgp(atot,gref,ptot)
1388 
1389  l = maxval(pgrid%xyz(xcoord,1:pgrid%nVert))
1390  ro = maxval(pgrid%xyz(ycoord,1:pgrid%nVert))
1391 
1392  im = max(pmixtinput%prepIntVal1,1)
1393  in = max(pmixtinput%prepIntVal2,1)
1394  iq = max(pmixtinput%prepIntVal3,1)
1395  ibc = max(min(pmixtinput%prepIntVal4,1),0)
1396 
1397  const = max(pmixtinput%prepRealVal1,0.0_rfreal)
1398 
1399  CALL rflu_jyzom(im,iq,dummyreal,etaqm,dummyreal,dummyreal)
1400 
1401  omega = atot*sqrt((in*global%pi/l)**2 + (etaqm/ro)**2)
1402 
1403  IF ( global%verbLevel > verbose_low ) THEN
1404  WRITE(stdout,'(A,5X,A,1X,I2)' ) solver_name, &
1405  'Boundary condition:',ibc
1406  WRITE(stdout,'(A,5X,A,3(1X,I2))') solver_name, &
1407  'Mode:',im,in,iq
1408  WRITE(stdout,'(A,5X,A,1X,E13.6)') solver_name, &
1409  'Total density (kg/m^3): ',dtot
1410  WRITE(stdout,'(A,5X,A,1X,E13.6)') solver_name, &
1411  'Total pressure (N/m^2): ',ptot
1412  WRITE(stdout,'(A,5X,A,1X,E13.6)') solver_name, &
1413  'Angular frequency (rad/s):',omega
1414  WRITE(stdout,'(A,5X,A,1X,E13.6)') solver_name, &
1415  'Constant (-): ',const
1416  END IF ! global%verbLevel
1417 
1418  DO icg = 1,pgrid%nCellsTot
1419  x = pgrid%cofg(xcoord,icg)
1420  y = pgrid%cofg(ycoord,icg)
1421  z = pgrid%cofg(zcoord,icg)
1422 
1423  CALL rflu_computeexactflowpacoust(global,z,y,x,global%currentTime, &
1424  l,ro,ibc,im,in,iq,etaqm,omega, &
1425  dtot,ptot,atot,const,d,u,v,w,p)
1426 
1427  pcv(cv_mixt_dens,icg) = d
1428  pcv(cv_mixt_xmom,icg) = d*u
1429  pcv(cv_mixt_ymom,icg) = d*v
1430  pcv(cv_mixt_zmom,icg) = d*w
1431  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,gref,p,u,v,w)
1432  END DO ! icg
1433 
1434 ! ------------------------------------------------------------------------------
1435 ! Ringleb flow. NOTE this problem is two-dimensional and assumed to lie in
1436 ! the x-y plane and that the exact solution is restricted to gamma = 1.4.
1437 ! ------------------------------------------------------------------------------
1438 
1439  CASE ( "ringleb" )
1440  CALL rflu_getparamshardcoderingleb(ptot,ttot)
1441 
1442  DO icg = 1,pgrid%nCellsTot
1443  x = pgrid%cofg(xcoord,icg)
1444  y = pgrid%cofg(ycoord,icg)
1445 
1446  CALL rflu_computeexactflowringleb(x,y,gcref,ptot,ttot,d,u,v,w,p)
1447 
1448  pcv(cv_mixt_dens,icg) = d
1449  pcv(cv_mixt_xmom,icg) = d*u
1450  pcv(cv_mixt_ymom,icg) = d*v
1451  pcv(cv_mixt_zmom,icg) = d*w
1452  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,1.4_rfreal,p,u,v,w)
1453  END DO ! icg
1454 
1455 ! ------------------------------------------------------------------------------
1456 ! Shock-bubble interaction (Quirk and Karni 1996)
1457 ! ------------------------------------------------------------------------------
1458 
1459  CASE ( "ShockBubble" )
1460  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
1461  CALL errorstop(global,err_gasmodel_invalid,__line__, &
1462  'Case initialization only valid with gas-liq model.')
1463  END IF ! pRegion%mixtInput%gasModel
1464 
1465  IF ( pregion%specInput%nSpecies /= 2 ) THEN
1466  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
1467  pregion%specInput%nSpecies
1468  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
1469  trim(errorstring))
1470  END IF ! pRegion%specInput%nSpecies
1471 
1472  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
1473  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
1474  rgas)
1475  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
1476  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
1477  rvap)
1478  cvl = global%refCvLiq
1479 
1480  DO icg = 1,pgrid%nCellsTot
1481  x = pgrid%cofg(xcoord,icg)
1482  y = pgrid%cofg(ycoord,icg)
1483 
1484  IF ( ((x-0.085_rfreal)**2 + y**2) <= 6.25e-04_rfreal )THEN
1485  rl = 1000.0_rfreal
1486  pl = 0.0_rfreal
1487  rg = 1.0_rfreal
1488  pg = 0.0_rfreal
1489  rv = 0.232127_rfreal
1490  pv = 1.0_rfreal
1491 
1492  d = rl*pl + rg*pg + rv*pv
1493  u = 0.0_rfreal
1494  v = 0.0_rfreal
1495  w = 0.0_rfreal
1496  t = 273.0_rfreal
1497  ELSE IF ( x < 0.050_rfreal ) THEN
1498  rl = 1000.0_rfreal
1499  pl = 0.0_rfreal
1500  rg = 1.75666_rfreal
1501  pg = 1.0_rfreal
1502  rv = 0.0_rfreal
1503  pv = 0.0_rfreal
1504 
1505  d = rl*pl + rg*pg + rv*pv
1506  u = 110.49_rfreal
1507  v = 0.0_rfreal
1508  w = 0.0_rfreal
1509  t = 311.366_rfreal
1510  ELSE
1511  rl = 1000_rfreal
1512  pl = 0.0_rfreal
1513  rg = 1.2763_rfreal
1514  pg = 1.0_rfreal
1515  rv = 0.0_rfreal
1516  pv = 0.0_rfreal
1517 
1518  d = rl*pl + rg*pg + rv*pv
1519  u = 0.0_rfreal
1520  v = 0.0_rfreal
1521  w = 0.0_rfreal
1522  t = 273.0_rfreal
1523  END IF ! x
1524 
1525  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
1526  vm2 = u*u + v*v + w*w
1527 
1528  pcv(cv_mixt_dens,icg) = d
1529  pcv(cv_mixt_xmom,icg) = d*u
1530  pcv(cv_mixt_ymom,icg) = d*v
1531  pcv(cv_mixt_zmom,icg) = d*w
1532  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
1533  END DO ! icg
1534 
1535 
1536 ! ------------------------------------------------------------------------------
1537 ! Skews diffracting shock
1538 ! ------------------------------------------------------------------------------
1539 
1540  CASE ( "skews_ms2p0","skews_ms3p0","skews_ms4p0" )
1541  DO icg = 1,pgrid%nCellsTot
1542  x = pgrid%cofg(xcoord,icg)
1543 
1544  IF ( x < pmixtinput%prepRealVal1 ) THEN
1545  d = pmixtinput%prepRealVal2
1546  u = pmixtinput%prepRealVal3
1547  v = 0.0_rfreal
1548  w = 0.0_rfreal
1549  p = pmixtinput%prepRealVal4
1550  ELSE
1551  d = pmixtinput%prepRealVal5
1552  u = 0.0_rfreal
1553  v = 0.0_rfreal
1554  w = 0.0_rfreal
1555  p = pmixtinput%prepRealVal6
1556  END IF ! z
1557 
1558  mw = pgv(gv_mixt_mol,indmol*icg)
1559  cp = pgv(gv_mixt_cp ,indcp *icg)
1560 
1561  gc = mixtperf_r_m(mw)
1562  g = mixtperf_g_cpr(cp,gc)
1563 
1564  pcv(cv_mixt_dens,icg) = d
1565  pcv(cv_mixt_xmom,icg) = d*u
1566  pcv(cv_mixt_ymom,icg) = d*v
1567  pcv(cv_mixt_zmom,icg) = d*w
1568  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1569  END DO ! icg
1570 
1571 ! ------------------------------------------------------------------------------
1572 ! Sommerfeld shock-particle-interaction
1573 ! ------------------------------------------------------------------------------
1574 
1575  CASE ( "somm_spi" )
1576  DO icg = 1,pgrid%nCellsTot
1577  x = pgrid%cofg(xcoord,icg)
1578 
1579  IF ( x < pmixtinput%prepRealVal1 ) THEN
1580  d = pmixtinput%prepRealVal2
1581  u = pmixtinput%prepRealVal3
1582  v = 0.0_rfreal
1583  w = 0.0_rfreal
1584  p = pmixtinput%prepRealVal4
1585  ELSE
1586  d = pmixtinput%prepRealVal5
1587  u = 0.0_rfreal
1588  v = 0.0_rfreal
1589  w = 0.0_rfreal
1590  p = pmixtinput%prepRealVal6
1591  END IF ! z
1592 
1593  mw = pgv(gv_mixt_mol,indmol*icg)
1594  cp = pgv(gv_mixt_cp ,indcp *icg)
1595 
1596  gc = mixtperf_r_m(mw)
1597  g = mixtperf_g_cpr(cp,gc)
1598 
1599  pcv(cv_mixt_dens,icg) = d
1600  pcv(cv_mixt_xmom,icg) = d*u
1601  pcv(cv_mixt_ymom,icg) = d*v
1602  pcv(cv_mixt_zmom,icg) = d*w
1603  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1604  END DO ! icg
1605 
1606 ! ------------------------------------------------------------------------------
1607 ! Sphere shock diffraction
1608 ! ------------------------------------------------------------------------------
1609 
1610  CASE ( "sphds" )
1611  DO icg = 1,pgrid%nCellsTot
1612  x = pgrid%cofg(xcoord,icg)
1613 
1614  IF ( x < pmixtinput%prepRealVal1 ) THEN
1615  d = pmixtinput%prepRealVal2
1616  u = pmixtinput%prepRealVal3
1617  v = 0.0_rfreal
1618  w = 0.0_rfreal
1619  p = pmixtinput%prepRealVal4
1620  ELSE
1621  d = pmixtinput%prepRealVal5
1622  u = 0.0_rfreal
1623  v = 0.0_rfreal
1624  w = 0.0_rfreal
1625  p = pmixtinput%prepRealVal6
1626  END IF ! x
1627 
1628  mw = pgv(gv_mixt_mol,indmol*icg)
1629  cp = pgv(gv_mixt_cp ,indcp *icg)
1630 
1631  gc = mixtperf_r_m(mw)
1632  g = mixtperf_g_cpr(cp,gc)
1633 
1634  pcv(cv_mixt_dens,icg) = d
1635  pcv(cv_mixt_xmom,icg) = d*u
1636  pcv(cv_mixt_ymom,icg) = d*v
1637  pcv(cv_mixt_zmom,icg) = d*w
1638  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1639  END DO ! icg
1640 
1641 ! ------------------------------------------------------------------------------
1642 ! Shocktubes
1643 ! ------------------------------------------------------------------------------
1644 
1645 ! ----- Generic 1d/2d ----------------------------------------------------------
1646 
1647  CASE ( "stg1d","stg2d" )
1648  DO icg = 1,pgrid%nCellsTot
1649  x = pgrid%cofg(xcoord,icg)
1650 
1651  IF ( x < pmixtinput%prepRealVal1 ) THEN
1652  d = pmixtinput%prepRealVal2
1653  u = pmixtinput%prepRealVal3
1654  v = 0.0_rfreal
1655  w = 0.0_rfreal
1656  p = pmixtinput%prepRealVal4
1657  ELSE IF ( (x >= pmixtinput%prepRealVal1 ) .AND. &
1658  (x < pmixtinput%prepRealVal11) ) THEN
1659  d = pmixtinput%prepRealVal5
1660  u = pmixtinput%prepRealVal6
1661  v = 0.0_rfreal
1662  w = 0.0_rfreal
1663  p = pmixtinput%prepRealVal7
1664  ELSE
1665  d = pmixtinput%prepRealVal12
1666  u = pmixtinput%prepRealVal13
1667  v = 0.0_rfreal
1668  w = 0.0_rfreal
1669  p = pmixtinput%prepRealVal14
1670  END IF ! x
1671 
1672  mw = pgv(gv_mixt_mol,indmol*icg)
1673  cp = pgv(gv_mixt_cp ,indcp *icg)
1674 
1675  gc = mixtperf_r_m(mw)
1676  g = mixtperf_g_cpr(cp,gc)
1677 
1678  pcv(cv_mixt_dens,icg) = d
1679  pcv(cv_mixt_xmom,icg) = d*u
1680  pcv(cv_mixt_ymom,icg) = d*v
1681  pcv(cv_mixt_zmom,icg) = d*w
1682  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1683  END DO ! icg
1684 
1685 ! ----- Sod case 1 -------------------------------------------------------------
1686 
1687  CASE ( "st_sod1","st_sod1_mp2" )
1688  DO icg = 1,pgrid%nCellsTot
1689  x = pgrid%cofg(xcoord,icg)
1690 
1691  IF ( x < 0.5_rfreal ) THEN
1692  d = 1.0_rfreal
1693  u = 0.0_rfreal
1694  v = 0.0_rfreal
1695  w = 0.0_rfreal
1696  p = 1.0e+5_rfreal
1697  ELSE
1698  d = 0.125_rfreal
1699  u = 0.0_rfreal
1700  v = 0.0_rfreal
1701  w = 0.0_rfreal
1702  p = 1.0e+4_rfreal
1703  END IF ! x
1704 
1705  mw = pgv(gv_mixt_mol,indmol*icg)
1706  cp = pgv(gv_mixt_cp ,indcp *icg)
1707 
1708  gc = mixtperf_r_m(mw)
1709  g = mixtperf_g_cpr(cp,gc)
1710 
1711  pcv(cv_mixt_dens,icg) = d
1712  pcv(cv_mixt_xmom,icg) = d*u
1713  pcv(cv_mixt_ymom,icg) = d*v
1714  pcv(cv_mixt_zmom,icg) = d*w
1715  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1716  END DO ! icg
1717 
1718 ! ----- Sod case 2 -------------------------------------------------------------
1719 
1720  CASE ( "st_sod2","st_sod2_mp2" )
1721  DO icg = 1,pgrid%nCellsTot
1722  x = pgrid%cofg(xcoord,icg)
1723 
1724  IF ( x < 0.0_rfreal ) THEN
1725  d = 1.0_rfreal
1726  u = 0.0_rfreal
1727  v = 0.0_rfreal
1728  w = 0.0_rfreal
1729  p = 1.0e+5_rfreal
1730  ELSE
1731  d = 0.01_rfreal
1732  u = 0.0_rfreal
1733  v = 0.0_rfreal
1734  w = 0.0_rfreal
1735  p = 1.0e+3_rfreal
1736  END IF ! x
1737 
1738  mw = pgv(gv_mixt_mol,indmol*icg)
1739  cp = pgv(gv_mixt_cp ,indcp *icg)
1740 
1741  gc = mixtperf_r_m(mw)
1742  g = mixtperf_g_cpr(cp,gc)
1743 
1744  pcv(cv_mixt_dens,icg) = d
1745  pcv(cv_mixt_xmom,icg) = d*u
1746  pcv(cv_mixt_ymom,icg) = d*v
1747  pcv(cv_mixt_zmom,icg) = d*w
1748  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1749  END DO ! icg
1750 
1751 ! ------------------------------------------------------------------------------
1752 ! Supersonic vortex flow. NOTE this problem is two-dimensional and assumed
1753 ! to lie in the x-y plane.
1754 ! ------------------------------------------------------------------------------
1755 
1756  CASE ( "ssvorth20x5l1" ,"ssvortp20x5l1" , &
1757  "ssvorth20x5l3" ,"ssvortp20x5l3" , &
1758  "ssvorth40x10l1" ,"ssvortp40x10l1" , &
1759  "ssvorth40x10l3" ,"ssvortp40x10l3" , &
1760  "ssvorth80x20l1" ,"ssvortp80x20l1" , &
1761  "ssvorth80x20l3" ,"ssvortp80x20l3" , &
1762  "ssvorth160x40l1" ,"ssvortp160x40l1" , &
1763  "ssvorth160x40l3" ,"ssvortp160x40l3" , &
1764  "ssvorth320x80l1" ,"ssvortp320x80l1" , &
1765  "ssvorth320x80l3" ,"ssvortp320x80l3" , &
1766  "ssvorth640x160l1","ssvortp640x160l1", &
1767  "ssvorth640x160l3","ssvortp640x160l3" )
1768  CALL rflu_getparamshardcodessvortex(ri,mi,ptot,ttot)
1769 
1770  DO icg = 1,pgrid%nCellsTot
1771  x = pgrid%cofg(xcoord,icg)
1772  y = pgrid%cofg(ycoord,icg)
1773 
1774  CALL rflu_computeexactflowssvortex(x,y,gref,gcref,ri,mi,ptot, &
1775  ttot,d,u,v,w,p)
1776 
1777  pcv(cv_mixt_dens,icg) = d
1778  pcv(cv_mixt_xmom,icg) = d*u
1779  pcv(cv_mixt_ymom,icg) = d*v
1780  pcv(cv_mixt_zmom,icg) = d*w
1781  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,gref,p,u,v,w)
1782  END DO ! icg
1783 
1784 ! ------------------------------------------------------------------------------
1785 ! Multiphase Riemann problem: Two rarefactions (Toro case 1)
1786 ! ------------------------------------------------------------------------------
1787 
1788  CASE ( "Two_Rarefaction" )
1789  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
1790  CALL errorstop(global,err_gasmodel_invalid,__line__, &
1791  'Case initialization only valid with gas-liq model.')
1792  END IF ! pRegion%mixtInput%gasModel
1793 
1794  IF ( pregion%specInput%nSpecies /= 2 ) THEN
1795  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
1796  pregion%specInput%nSpecies
1797  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
1798  trim(errorstring))
1799  END IF ! pRegion%specInput%nSpecies
1800 
1801  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
1802  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
1803  rgas)
1804  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
1805  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
1806  rvap)
1807  cvl = global%refCvLiq
1808 
1809  DO icg = 1,pgrid%nCellsTot
1810  x = pgrid%cofg(xcoord,icg)
1811 
1812  IF ( x < 0.5_rfreal ) THEN
1813  rl = 1000_rfreal
1814  pl = 0.1_rfreal
1815  rg = 1.2342_rfreal
1816  pg = 0.9_rfreal
1817  rv = 0.0_rfreal
1818  pv = 0.0_rfreal
1819 
1820  d = rl*pl + rg*pg + rv*pv
1821  u = -200.0_rfreal
1822  v = 0.0_rfreal
1823  w = 0.0_rfreal
1824  t = 273.0_rfreal
1825  ELSE
1826  rl = 1000_rfreal
1827  pl = 0.1_rfreal
1828  rg = 1.2342_rfreal
1829  pg = 0.9_rfreal
1830  rv = 0.0_rfreal
1831  pv = 0.0_rfreal
1832 
1833  d = rl*pl + rg*pg + rv*pv
1834  u = 200.0_rfreal
1835  v = 0.0_rfreal
1836  w = 0.0_rfreal
1837  t = 273.0_rfreal
1838  END IF ! x
1839 
1840  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
1841  vm2 = u*u+v*v+w*w
1842 
1843  pcv(cv_mixt_dens,icg) = d
1844  pcv(cv_mixt_xmom,icg) = d*u
1845  pcv(cv_mixt_ymom,icg) = d*v
1846  pcv(cv_mixt_zmom,icg) = d*w
1847  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
1848  END DO ! icg
1849 
1850 ! ------------------------------------------------------------------------------
1851 ! Simple volcano model
1852 ! ------------------------------------------------------------------------------
1853 
1854  CASE ( "volcmod2dv3" )
1855  DO icg = 1,pgrid%nCellsTot
1856  y = pgrid%cofg(ycoord,icg)
1857 
1858  IF ( y < -500.0_rfreal ) THEN
1859  u = 0.0_rfreal
1860  v = 0.0_rfreal
1861  w = 0.0_rfreal
1862  p = 12.5e+6_rfreal
1863  t = 600.0_rfreal
1864  ELSE
1865  u = 0.0_rfreal
1866  v = 0.0_rfreal
1867  w = 0.0_rfreal
1868  p = 0.87e+5_rfreal
1869  t = 288.15_rfreal
1870  END IF ! y
1871 
1872  mw = pgv(gv_mixt_mol,indmol*icg)
1873  cp = pgv(gv_mixt_cp ,indcp *icg)
1874 
1875  gc = mixtperf_r_m(mw)
1876  g = mixtperf_g_cpr(cp,gc)
1877 
1878  d = mixtperf_d_prt(p,gc,t)
1879 
1880  pcv(cv_mixt_dens,icg) = d
1881  pcv(cv_mixt_xmom,icg) = d*u
1882  pcv(cv_mixt_ymom,icg) = d*v
1883  pcv(cv_mixt_zmom,icg) = d*w
1884  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1885  END DO ! icg
1886 
1887 ! ------------------------------------------------------------------------------
1888 ! Vortex test
1889 ! ------------------------------------------------------------------------------
1890 
1891  CASE ( "vort" )
1892 
1893  rc_l = 0.15_rfreal
1894  a_cl = -0.0005_rfreal
1895  uo_c = 1.1_rfreal
1896 ! uo_c = 0.3_RFREAL
1897 
1898  rc = rc_l*pmixtinput%prepRealVal1
1899  c = pmixtinput%prepRealVal2
1900  a = a_cl*c*pmixtinput%prepRealVal1
1901  uo = uo_c*c
1902 
1903  ro = pmixtinput%prepRealVal3
1904  po = pmixtinput%prepRealVal4
1905 
1906  DO icg = 1,pgrid%nCellsTot
1907  x = pgrid%cofg(xcoord,icg)
1908  y = pgrid%cofg(ycoord,icg)
1909 
1910  radius = sqrt(x*x+y*y)
1911 
1912  mw = pgv(gv_mixt_mol,indmol*icg)
1913  cp = pgv(gv_mixt_cp ,indcp *icg)
1914 
1915  gc = mixtperf_r_m(mw)
1916  g = mixtperf_g_cpr(cp,gc)
1917 
1918  IF ( radius <= rc ) THEN
1919  psi = a*exp(-(x*x+y*y)/(2.0_rfreal*rc*rc))
1920  ELSE
1921  psi = 0.0_rfreal
1922  END IF ! radius
1923 
1924  d = ro
1925  u = uo - y*psi/(d*rc*rc)
1926  v = x*psi/(d*rc*rc)
1927  w = 0.0_rfreal
1928  p = po + (d*a*psi)/(rc*rc)
1929 
1930  pcv(cv_mixt_dens,icg) = d
1931  pcv(cv_mixt_xmom,icg) = d*u
1932  pcv(cv_mixt_ymom,icg) = d*v
1933  pcv(cv_mixt_zmom,icg) = d*w
1934  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1935  END DO ! icg
1936 
1937 ! ------------------------------------------------------------------------------
1938 ! Solid-body vortex
1939 ! ------------------------------------------------------------------------------
1940 
1941  CASE ( "vortex_part" )
1942  DO icg = 1,pgrid%nCellsTot
1943  x = pgrid%cofg(xcoord,icg)
1944  y = pgrid%cofg(ycoord,icg)
1945 
1946  mw = pgv(gv_mixt_mol,indmol*icg)
1947  cp = pgv(gv_mixt_cp ,indcp *icg)
1948 
1949  gc = mixtperf_r_m(mw)
1950  g = mixtperf_g_cpr(cp,gc)
1951 
1952  d = 1.0_rfreal
1953  u = y
1954  v = -x
1955  w = 0.0_rfreal
1956  p = 1.0e+5_rfreal
1957 
1958  pcv(cv_mixt_dens,icg) = d
1959  pcv(cv_mixt_xmom,icg) = d*u
1960  pcv(cv_mixt_ymom,icg) = d*v
1961  pcv(cv_mixt_zmom,icg) = d*w
1962  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1963  END DO ! icg
1964 
1965 ! ------------------------------------------------------------------------------
1966 ! Woodward-Colella ramp. NOTE the initial condition is one-dimensional
1967 ! and assumed to lie along the x-axis.
1968 ! ------------------------------------------------------------------------------
1969 
1970  CASE ( "wcramp","wcrampsc","wcrampsm" )
1971  DO icg = 1,pgrid%nCellsTot
1972  x = pgrid%cofg(xcoord,icg)
1973 
1974  IF ( x < -0.1_rfreal ) THEN
1975  d = 8.0_rfreal
1976  u = 2608.87906904_rfreal
1977  v = 0.0_rfreal
1978  w = 0.0_rfreal
1979  p = 116.5e+5_rfreal
1980  ELSE
1981  d = 1.4_rfreal
1982  u = 0.0_rfreal
1983  v = 0.0_rfreal
1984  w = 0.0_rfreal
1985  p = 1.0e+5_rfreal
1986  END IF ! x
1987 
1988  mw = pgv(gv_mixt_mol,indmol*icg)
1989  cp = pgv(gv_mixt_cp ,indcp *icg)
1990 
1991  gc = mixtperf_r_m(mw)
1992  g = mixtperf_g_cpr(cp,gc)
1993 
1994  pcv(cv_mixt_dens,icg) = d
1995  pcv(cv_mixt_xmom,icg) = d*u
1996  pcv(cv_mixt_ymom,icg) = d*v
1997  pcv(cv_mixt_zmom,icg) = d*w
1998  pcv(cv_mixt_ener,icg) = d*mixtperf_eo_dgpuvw(d,g,p,u,v,w)
1999  END DO ! icg
2000 
2001 ! ------------------------------------------------------------------------------
2002 ! Gas-liquid mixture: Sod shocktube
2003 ! ------------------------------------------------------------------------------
2004 
2005  CASE ( "2DShock001" )
2006  IF ( pregion%mixtInput%gasModel /= gas_model_mixt_gasliq ) THEN
2007  CALL errorstop(global,err_gasmodel_invalid,__line__, &
2008  'Case initialization only valid with gas-liq model.')
2009  END IF ! pRegion%mixtInput%gasModel
2010 
2011  IF ( pregion%specInput%nSpecies /= 2 ) THEN
2012  WRITE(errorstring,'(A,1X,I2)') 'Should be:', &
2013  pregion%specInput%nSpecies
2014  CALL errorstop(global,err_spec_nspec_invalid,__line__, &
2015  trim(errorstring))
2016  END IF ! pRegion%specInput%nSpecies
2017 
2018  rgas = mixtperf_r_m(pregion%specInput%specType(1)%pMaterial%molw)
2019  cvg = mixtperf_cv_cpr(pregion%specInput%specType(1)%pMaterial%spht, &
2020  rgas)
2021  rvap = mixtperf_r_m(pregion%specInput%specType(2)%pMaterial%molw)
2022  cvv = mixtperf_cv_cpr(pregion%specInput%specType(2)%pMaterial%spht, &
2023  rvap)
2024  cvl = global%refCvLiq
2025 
2026  DO icg = 1,pgrid%nCellsTot
2027  x = pgrid%cofg(xcoord,icg)
2028 
2029  IF ( x < 0.5_rfreal ) THEN
2030  rl = 0.0_rfreal
2031  pl = 0.0_rfreal
2032  rg = 1.0_rfreal
2033  pg = 1.0_rfreal
2034  rv = 0.0_rfreal
2035  pv = 0.0_rfreal
2036 
2037  d = rl*pl + rg*pg + rv*pv
2038  u = 0.0_rfreal
2039  v = 0.0_rfreal
2040  w = 0.0_rfreal
2041  t = 348.43_rfreal
2042  ELSE
2043  rl = 0.0_rfreal
2044  pl = 0.0_rfreal
2045  rg = 0.125_rfreal
2046  pg = 1.0_rfreal
2047  rv = 0.0_rfreal
2048  pv = 0.0_rfreal
2049 
2050  d = rl*pl + rg*pg + rv*pv
2051  u = 0.0_rfreal
2052  v = 0.0_rfreal
2053  w = 0.0_rfreal
2054  t = 278.7_rfreal
2055  END IF ! x
2056 
2057  cvm = (rl*pl*cvl + rg*pg*cvg + rv*pv*cvv)/d
2058  vm2 = u*u+v*v+w*w
2059 
2060  pcv(cv_mixt_dens,icg) = d
2061  pcv(cv_mixt_xmom,icg) = d*u
2062  pcv(cv_mixt_ymom,icg) = d*v
2063  pcv(cv_mixt_zmom,icg) = d*w
2064  pcv(cv_mixt_ener,icg) = d*mixtgasliq_eo_cvmtvm2(cvm,t,vm2)
2065  END DO ! icg
2066 
2067 ! ------------------------------------------------------------------------------
2068 ! Default - must be due to input error
2069 ! ------------------------------------------------------------------------------
2070 
2071  CASE default
2072  CALL errorstop(global,err_reached_default,__line__)
2073  END SELECT ! global%casename
2074 
2075 ! ==============================================================================
2076 ! Default
2077 ! ==============================================================================
2078 
2079  CASE default
2080  CALL errorstop(global,err_reached_default,__line__)
2081  END SELECT ! pMixtInput%fluidModel
2082 
2083 ! ******************************************************************************
2084 ! End
2085 ! ******************************************************************************
2086 
2087  IF ( global%verbLevel > verbose_none ) THEN
2088  WRITE(stdout,'(A,1X,A)') solver_name, &
2089  'Initializing flow field from hard code done.'
2090  END IF ! global%verbLevel
2091 
2092  CALL deregisterfunction(global)
2093 
2094 END SUBROUTINE rflu_initflowhardcode
2095 
2096 ! ******************************************************************************
2097 !
2098 ! RCS Revision history:
2099 !
2100 ! $Log: RFLU_InitFlowHardCode.F90,v $
2101 ! Revision 1.38 2008/12/06 08:44:56 mtcampbe
2102 ! Updated license.
2103 !
2104 ! Revision 1.37 2008/11/19 22:18:06 mtcampbe
2105 ! Added Illinois Open Source License/Copyright
2106 !
2107 ! Revision 1.36 2007/04/05 01:13:16 haselbac
2108 ! Added stg1d, modified code to allow 2nd interface
2109 !
2110 ! Revision 1.35 2007/02/27 13:15:00 haselbac
2111 ! Added init for stg1d
2112 !
2113 ! Revision 1.34 2007/02/16 20:01:05 haselbac
2114 ! Added init for somm_spi case
2115 !
2116 ! Revision 1.33 2006/08/19 19:44:09 haselbac
2117 ! Significant clean-up and cosmetic changes
2118 !
2119 ! Revision 1.32 2006/08/19 15:41:11 mparmar
2120 ! Added initializations for vortexNSCBC, farf, nscbc[1-8], vort
2121 !
2122 ! Revision 1.31 2006/05/06 16:50:56 haselbac
2123 ! Added MP versions of kjet2 cases
2124 !
2125 ! Revision 1.30 2006/04/20 20:57:35 haselbac
2126 ! Added kjet2v5 case
2127 !
2128 ! Revision 1.29 2006/04/19 19:26:23 haselbac
2129 ! Fixed order of CASE statements
2130 !
2131 ! Revision 1.28 2006/04/13 18:10:12 haselbac
2132 ! Added treatment of Culick flow
2133 !
2134 ! Revision 1.27 2006/03/30 20:53:04 haselbac
2135 ! Changed ShockBubble hard-code, cosmetics
2136 !
2137 ! Revision 1.26 2006/03/28 00:36:36 haselbac
2138 ! Added kjet2v4 case
2139 !
2140 ! Revision 1.25 2006/03/26 20:22:21 haselbac
2141 ! Added cases for GL model
2142 !
2143 ! Revision 1.24 2006/03/08 23:39:19 haselbac
2144 ! Added kjet2v3 case
2145 !
2146 ! Revision 1.23 2006/01/06 22:16:49 haselbac
2147 ! Added routines to init linear and trig fn for grad testing
2148 !
2149 ! Revision 1.22 2005/12/07 17:58:19 fnajjar
2150 ! Bug fix for gas constant defs in vortex_part
2151 !
2152 ! Revision 1.21 2005/12/07 16:58:37 haselbac
2153 ! Added init for vortex_part
2154 !
2155 ! Revision 1.20 2005/11/11 16:58:08 haselbac
2156 ! Added init for generic 2d shocktube
2157 !
2158 ! Revision 1.19 2005/11/10 16:58:33 haselbac
2159 ! Bug fix for gmpjet
2160 !
2161 ! Revision 1.18 2005/11/10 02:42:20 haselbac
2162 ! Added support for variable props, new cases
2163 !
2164 ! Revision 1.17 2005/10/09 15:12:17 haselbac
2165 ! Added 2d C0 case
2166 !
2167 ! Revision 1.16 2005/09/28 22:54:26 mparmar
2168 ! Changed cylds init to use prepRealValx
2169 !
2170 ! Revision 1.15 2005/09/13 21:38:15 haselbac
2171 ! Changed sphds init to use prepRealValx
2172 !
2173 ! Revision 1.14 2005/08/25 16:22:36 haselbac
2174 ! Bug fix in ds v6 init
2175 !
2176 ! Revision 1.13 2005/08/25 03:42:48 haselbac
2177 ! Added v6 to ds cases
2178 !
2179 ! Revision 1.12 2005/08/24 01:37:04 haselbac
2180 ! Added v5 for ds cases
2181 !
2182 ! Revision 1.11 2005/08/21 16:03:06 haselbac
2183 ! Added v4 to ds cases
2184 !
2185 ! Revision 1.10 2005/08/17 23:00:00 haselbac
2186 ! Added v3 versions of diffracting shock cases
2187 !
2188 ! Revision 1.9 2005/08/08 00:12:46 haselbac
2189 ! Modified ds case initialization, now use prepRealValx
2190 !
2191 ! Revision 1.8 2005/08/06 00:56:42 haselbac
2192 ! Added v2 to ds case names
2193 !
2194 ! Revision 1.7 2005/07/05 19:47:58 mparmar
2195 ! Added init for diffracting shock over sphere
2196 !
2197 ! Revision 1.6 2005/07/04 14:58:37 haselbac
2198 ! Added init for mvsh case
2199 !
2200 ! Revision 1.5 2005/07/01 16:21:36 haselbac
2201 ! Added init for kjet2
2202 !
2203 ! Revision 1.4 2005/06/14 00:57:45 haselbac
2204 ! Changed init for Skews case
2205 !
2206 ! Revision 1.3 2005/04/29 12:51:24 haselbac
2207 ! Adapted to changes in interface of RFLU_ComputeExactFlowPAcoust
2208 !
2209 ! Revision 1.2 2005/04/20 14:45:37 haselbac
2210 ! Extended init of pipeacoust case
2211 !
2212 ! Revision 1.1 2005/04/15 15:08:16 haselbac
2213 ! Initial revision
2214 !
2215 ! ******************************************************************************
2216 
2217 
2218 
2219 
2220 
2221 
2222 
subroutine, public rflu_getparamshardcodessvortex(ri, Mi, pTot, tTot)
subroutine, public rflu_computeexactflowculick(global, x, y, z, d, u, v, w, p)
unsigned char r() const
Definition: Color.h:68
real(rfreal) function mixtperf_r_m(M)
Definition: MixtPerf_R.F90:54
double ymin() const
const NT & d
void int int REAL REAL * y
Definition: read.cpp:74
subroutine, public rflu_computeexactflowpacoust(global, x, y, z, t, L, ro, iBc, im, in, iq, etaqm, omega, dTot, pTot, aTot, const, d, u, v, w, p)
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
double xmin() const
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
double zmin() const
real(rfreal) function mixtperf_d_prt(P, R, T)
Definition: MixtPerf_D.F90:71
double sqrt(double d)
Definition: double.h:73
subroutine, public rflu_jyzom(N, M, RJ0M, RJ1M, RY0M, RY1M)
real(rfreal) function mixtperf_r_cpg(Cp, G)
Definition: MixtPerf_R.F90:39
RT c() const
Definition: Line_2.h:150
*********************************************************************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
void int int int REAL REAL REAL * z
Definition: write.cpp:76
NT & sin
void int int REAL * x
Definition: read.cpp:74
real(rfreal) function mixtperf_d_cgp(C, G, P)
Definition: MixtPerf_D.F90:40
subroutine, public rflu_getparamshardcoderingleb(pTot, tTot)
subroutine, public rflu_computeexactflowringleb(x, y, rGas, pTot, tTot, d, u, v, w, p)
subroutine, public rflu_computeexactflowssvortex(x, y, gGas, rGas, ri, Mi, pTot, tTot, d, u, v, w, p)
real(rfreal) function mixtgasliq_eo_cvmtvm2(Cvm, T, Vm2)
Vector_n min(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:346
subroutine rflu_initflowhardcode(pRegion)
subroutine, public rflu_getparamshardcodeproudman(dInc, mInj, vInj, pTot)
subroutine, public rflu_setexactflowtrig(global, nx, ny, nz, x, y, z, iVar, var, gx, gy, gz)
real(rfreal) function mixtperf_eo_dgpuvw(D, G, P, U, V, W)
Definition: MixtPerf_E.F90:40
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_getparamshardcodepacoust(pTot, aTot)
real(rfreal) function mixtperf_g_cpr(Cp, R)
Definition: MixtPerf_G.F90:39
subroutine, public rflu_setexactflowlinear(x, y, z, iVar, var, gx, gy, gz)
NT & cos
real(rfreal) function mixtperf_cv_cpr(Cp, R)
Definition: MixtPerf_Cv.F90:39
CImg< T > & atan2(const CImg< t > &img)
Compute the arc-tangent of each pixel.
Definition: CImg.h:12671
real(rfreal) function mixtperf_p_drt(D, R, T)
Definition: MixtPerf_P.F90:54
RT a() const
Definition: Line_2.h:140
unsigned char g() const
Definition: Color.h:69
subroutine, public rflu_computeexactflowproudman(global, x, y, height, dInc, vInj, pTot, d, u, v, w, p)