@@ -29,6 +29,7 @@ module m_derived_variables
2929 s_derive_flux_limiter, &
3030 s_derive_vorticity_component, &
3131 s_derive_qm, &
32+ s_derive_liutex, &
3233 s_derive_numerical_schlieren_function, &
3334 s_compute_speed_of_sound, &
3435 s_finalize_derived_variables_module
@@ -80,21 +81,21 @@ contains
8081 ! s_compute_finite_difference_coefficients.
8182
8283 ! Allocating centered finite- difference coefficients in x- direction
83- if (omega_wrt(2 ) .or. omega_wrt(3 ) .or. schlieren_wrt) then
84+ if (omega_wrt(2 ) .or. omega_wrt(3 ) .or. schlieren_wrt .or. liutex_wrt ) then
8485 allocate (fd_coeff_x(- fd_number:fd_number, &
8586 - offset_x%beg:m + offset_x%end))
8687 end if
8788
8889 ! Allocating centered finite- difference coefficients in y- direction
89- if (omega_wrt(1 ) .or. omega_wrt(3 ) &
90+ if (omega_wrt(1 ) .or. omega_wrt(3 ) .or. liutex_wrt &
9091 .or. &
9192 (n > 0 .and. schlieren_wrt)) then
9293 allocate (fd_coeff_y(- fd_number:fd_number, &
9394 - offset_y%beg:n + offset_y%end))
9495 end if
9596
9697 ! Allocating centered finite- difference coefficients in z- direction
97- if (omega_wrt(1 ) .or. omega_wrt(2 ) &
98+ if (omega_wrt(1 ) .or. omega_wrt(2 ) .or. liutex_wrt &
9899 .or. &
99100 (p > 0 .and. schlieren_wrt)) then
100101 allocate (fd_coeff_z(- fd_number:fd_number, &
@@ -556,6 +557,135 @@ contains
556557
557558 end subroutine s_derive_qm
558559
560+ !> This subroutine gets as inputs the primitive variables. From those
561+ !! inputs, it proceeds to calculate the Liutex vector and its
562+ !! magnitude based on Xu et al. (2019).
563+ !! @param q_prim_vf Primitive variables
564+ !! @param liutex_mag Liutex magnitude
565+ !! @param liutex_axis Liutex axis
566+ impure subroutine s_derive_liutex(q_prim_vf, liutex_mag, liutex_axis)
567+ integer, parameter :: nm = 3
568+ type(scalar_field), &
569+ dimension(sys_size), &
570+ intent(in) :: q_prim_vf
571+
572+ real(wp), &
573+ dimension(-offset_x%beg:m + offset_x%end, &
574+ -offset_y%beg:n + offset_y%end, &
575+ -offset_z%beg:p + offset_z%end), &
576+ intent(out) :: liutex_mag !< Liutex magnitude
577+
578+ real(wp), &
579+ dimension(-offset_x%beg:m + offset_x%end, &
580+ -offset_y%beg:n + offset_y%end, &
581+ -offset_z%beg:p + offset_z%end, nm), &
582+ intent(out) :: liutex_axis !< Liutex rigid rotation axis
583+
584+ character, parameter :: ivl = ' N' !< compute left eigenvectors
585+ character, parameter :: ivr = ' V' !< compute right eigenvectors
586+ real(wp), dimension(nm, nm) :: vgt !< velocity gradient tensor
587+ real(wp), dimension(nm) :: lr, li !< real and imaginary parts of eigenvalues
588+ real(wp), dimension(nm, nm) :: vl, vr !< left and right eigenvectors
589+ integer, parameter :: lwork = 4*nm !< size of work array (4*nm recommended)
590+ real(wp), dimension(lwork) :: work !< work array
591+ integer :: info
592+
593+ real(wp), dimension(nm) :: eigvec !< real eigenvector
594+ real(wp) :: eigvec_mag !< magnitude of real eigenvector
595+ real(wp) :: omega_proj !< projection of vorticity on real eigenvector
596+ real(wp) :: lci !< imaginary part of complex eigenvalue
597+ real(wp) :: alpha
598+
599+ integer :: j, k, l, r, i !< Generic loop iterators
600+ integer :: idx
601+
602+ do l = -offset_z%beg, p + offset_z%end
603+ do k = -offset_y%beg, n + offset_y%end
604+ do j = -offset_x%beg, m + offset_x%end
605+
606+ ! Get velocity gradient tensor (VGT)
607+ vgt(:, :) = 0._wp
608+
609+ do r = -fd_number, fd_number
610+ do i = 1, 3
611+ ! d()/dx
612+ vgt(i, 1) = &
613+ vgt(i, 1) + &
614+ fd_coeff_x(r, j)* &
615+ q_prim_vf(mom_idx%beg + i - 1)%sf(r + j, k, l)
616+ ! d()/dy
617+ vgt(i, 2) = &
618+ vgt(i, 2) + &
619+ fd_coeff_y(r, k)* &
620+ q_prim_vf(mom_idx%beg + i - 1)%sf(j, r + k, l)
621+ ! d()/dz
622+ vgt(i, 3) = &
623+ vgt(i, 3) + &
624+ fd_coeff_z(r, l)* &
625+ q_prim_vf(mom_idx%beg + i - 1)%sf(j, k, r + l)
626+ end do
627+ end do
628+
629+ ! Call appropriate LAPACK routine based on precision
630+ #ifdef MFC_SINGLE_PRECISION
631+ call sgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
632+ #else
633+ call dgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
634+ #endif
635+
636+ ! Find real eigenvector
637+ idx = 1
638+ do r = 2, 3
639+ if (abs(li(r)) < abs(li(idx))) then
640+ idx = r
641+ end if
642+ end do
643+ eigvec = vr(:, idx)
644+
645+ ! Normalize real eigenvector if it is effectively non-zero
646+ eigvec_mag = sqrt(eigvec(1)**2._wp &
647+ + eigvec(2)**2._wp &
648+ + eigvec(3)**2._wp)
649+ if (eigvec_mag > sgm_eps) then
650+ eigvec = eigvec/eigvec_mag
651+ else
652+ eigvec = 0._wp
653+ end if
654+
655+ ! Compute vorticity projected on the eigenvector
656+ omega_proj = (vgt(3, 2) - vgt(2, 3))*eigvec(1) &
657+ + (vgt(1, 3) - vgt(3, 1))*eigvec(2) &
658+ + (vgt(2, 1) - vgt(1, 2))*eigvec(3)
659+
660+ ! As eigenvector can have +/- signs, we can choose the sign
661+ ! so that omega_proj is positive
662+ if (omega_proj < 0._wp) then
663+ eigvec = -eigvec
664+ omega_proj = -omega_proj
665+ end if
666+
667+ ! Find imaginary part of complex eigenvalue
668+ lci = li(mod(idx, 3) + 1)
669+
670+ ! Compute Liutex magnitude
671+ alpha = omega_proj**2._wp - 4._wp*lci**2._wp ! (2*alpha)^2
672+ if (alpha > 0._wp) then
673+ liutex_mag(j, k, l) = omega_proj - sqrt(alpha)
674+ else
675+ liutex_mag(j, k, l) = omega_proj
676+ end if
677+
678+ ! Compute Liutex axis
679+ liutex_axis(j, k, l, 1) = eigvec(1)
680+ liutex_axis(j, k, l, 2) = eigvec(2)
681+ liutex_axis(j, k, l, 3) = eigvec(3)
682+
683+ end do
684+ end do
685+ end do
686+
687+ end subroutine s_derive_liutex
688+
559689 !> This subroutine gets as inputs the conservative variables
560690 !! and density. From those inputs, it proceeds to calculate
561691 !! the values of the numerical Schlieren function, which are
0 commit comments