@@ -238,6 +238,11 @@ subroutine s_assign_patch_species_primitive_variables_bubbles(patch_id, j, k, l,
238238 real (kind (0d0 )) :: orig_pi_inf
239239 real (kind (0d0 )) :: muR, muV
240240
241+ real (kind (0d0 )), dimension (int (E_idx - mom_idx% beg)) :: vel ! < velocity
242+ real (kind (0d0 )) :: pres ! < pressure
243+ real (kind (0d0 )) :: x_centroid, y_centroid
244+ real (kind (0d0 )) :: epsilon, beta
245+
241246 real (kind (0d0 )), dimension (sys_size) :: orig_prim_vf ! <
242247 ! ! Vector to hold original values of cell for smoothing purposes
243248
@@ -557,6 +562,47 @@ subroutine s_assign_patch_species_primitive_variables_bubbles(patch_id, j, k, l,
557562 end if
558563 end do
559564 end if
565+
566+ if (patch_icpp(patch_id)% geometry == 6 ) then
567+ x_centroid = patch_icpp(patch_id)% x_centroid
568+ y_centroid = patch_icpp(patch_id)% y_centroid
569+ epsilon = patch_icpp(patch_id)% epsilon
570+ beta = patch_icpp(patch_id)% beta
571+
572+ ! Reference density, velocity, pressure and specific heat ratio
573+ ! function of the isentropic vortex patch
574+ rho = patch_icpp(patch_id)% rho
575+ vel = patch_icpp(patch_id)% vel
576+ pres = patch_icpp(patch_id)% pres
577+ gamma = patch_icpp(patch_id)% gamma
578+
579+ ! Density
580+ q_prim_vf(1 )% sf(j, k, 0 ) = &
581+ rho* (1d0 - (rho/ pres)* (epsilon/ (2d0 * pi))* &
582+ (epsilon/ (8d0 * beta* (gamma + 1d0 )* pi))* &
583+ exp (2d0 * beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
584+ - (y_cc(k) - y_centroid)** 2 )) &
585+ )** gamma
586+
587+ ! Velocity
588+ q_prim_vf(2 )% sf(j, k, 0 ) = &
589+ vel(1 ) - (y_cc(k) - y_centroid)* (epsilon/ (2d0 * pi))* &
590+ exp (beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
591+ - (y_cc(k) - y_centroid)** 2 ))
592+ q_prim_vf(3 )% sf(j, k, 0 ) = &
593+ vel(2 ) + (x_cc(j) - x_centroid)* (epsilon/ (2d0 * pi))* &
594+ exp (beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
595+ - (y_cc(k) - y_centroid)** 2 ))
596+
597+ ! Pressure
598+ q_prim_vf(4 )% sf(j, k, 0 ) = &
599+ pres* (1d0 - (rho/ pres)* (epsilon/ (2d0 * pi))* &
600+ (epsilon/ (8d0 * beta* (gamma + 1d0 )* pi))* &
601+ exp (2d0 * beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
602+ - (y_cc(k) - y_centroid)** 2 )) &
603+ )** (gamma + 1d0 )
604+
605+ end if
560606
561607 ! Updating the patch identities bookkeeping variable
562608 if (1d0 - eta < 1d-16 ) patch_id_fp(j, k, l) = patch_id
@@ -585,6 +631,7 @@ subroutine s_assign_patch_species_primitive_variables(patch_id, j, k, l, &
585631 integer , intent (INOUT ), dimension (0 :m, 0 :n, 0 :p) :: patch_id_fp
586632 type (scalar_field), dimension (1 :sys_size) :: q_prim_vf
587633 real (kind (0d0 )) :: eta ! <
634+
588635 integer , intent (IN ) :: j, k, l
589636
590637 real (kind (0d0 )) :: rho
@@ -597,6 +644,11 @@ subroutine s_assign_patch_species_primitive_variables(patch_id, j, k, l, &
597644 ! ! function, respectively, obtained from the combination of primitive
598645 ! ! variables of the current and smoothing patches
599646
647+ real (kind (0d0 )), dimension (int (E_idx - mom_idx% beg)) :: vel ! < velocity
648+ real (kind (0d0 )) :: pres ! < pressure
649+ real (kind (0d0 )) :: x_centroid, y_centroid
650+ real (kind (0d0 )) :: epsilon, beta
651+
600652 real (kind (0d0 )), dimension (sys_size) :: orig_prim_vf ! <
601653 ! Vector to hold original values of cell for smoothing purposes
602654
@@ -713,6 +765,47 @@ subroutine s_assign_patch_species_primitive_variables(patch_id, j, k, l, &
713765 q_prim_vf(i)% sf(j, k, l) = q_prim_vf(E_idx)% sf(j, k, l)
714766 end do
715767 end if
768+
769+ if (patch_icpp(patch_id)% geometry == 6 ) then
770+ x_centroid = patch_icpp(patch_id)% x_centroid
771+ y_centroid = patch_icpp(patch_id)% y_centroid
772+ epsilon = patch_icpp(patch_id)% epsilon
773+ beta = patch_icpp(patch_id)% beta
774+
775+ ! Reference density, velocity, pressure and specific heat ratio
776+ ! function of the isentropic vortex patch
777+ rho = patch_icpp(patch_id)% rho
778+ vel = patch_icpp(patch_id)% vel
779+ pres = patch_icpp(patch_id)% pres
780+ gamma = patch_icpp(patch_id)% gamma
781+
782+ ! Density
783+ q_prim_vf(1 )% sf(j, k, 0 ) = &
784+ rho* (1d0 - (rho/ pres)* (epsilon/ (2d0 * pi))* &
785+ (epsilon/ (8d0 * beta* (gamma + 1d0 )* pi))* &
786+ exp (2d0 * beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
787+ - (y_cc(k) - y_centroid)** 2 )) &
788+ )** gamma
789+
790+ ! Velocity
791+ q_prim_vf(2 )% sf(j, k, 0 ) = &
792+ vel(1 ) - (y_cc(k) - y_centroid)* (epsilon/ (2d0 * pi))* &
793+ exp (beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
794+ - (y_cc(k) - y_centroid)** 2 ))
795+ q_prim_vf(3 )% sf(j, k, 0 ) = &
796+ vel(2 ) + (x_cc(j) - x_centroid)* (epsilon/ (2d0 * pi))* &
797+ exp (beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
798+ - (y_cc(k) - y_centroid)** 2 ))
799+
800+ ! Pressure
801+ q_prim_vf(4 )% sf(j, k, 0 ) = &
802+ pres* (1d0 - (rho/ pres)* (epsilon/ (2d0 * pi))* &
803+ (epsilon/ (8d0 * beta* (gamma + 1d0 )* pi))* &
804+ exp (2d0 * beta* (1d0 - (x_cc(j) - x_centroid)** 2 &
805+ - (y_cc(k) - y_centroid)** 2 )) &
806+ )** (gamma + 1d0 )
807+
808+ end if
716809
717810 ! Updating the patch identities bookkeeping variable
718811 if (1d0 - eta < 1d-16 ) patch_id_fp(j, k, l) = patch_id
0 commit comments