curl calculation overestimated RMS due to factor 2 instead of one for DC component and Nyquist component
This commit is contained in:
parent
e88cedc6ae
commit
8f32d03a9e
|
@ -740,7 +740,7 @@ real(pReal) function utilities_curlRMS()
|
|||
-tensorField_fourier(l,1,i,j,k)*xi(2,i,j,k))*TWOPIIMG
|
||||
enddo
|
||||
utilities_curlRMS = utilities_curlRMS + &
|
||||
2.0_pReal*sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)
|
||||
2.0_pReal*sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)! Has somewhere a conj. complex counterpart. Therefore count it twice.
|
||||
enddo
|
||||
do l = 1_pInt, 3_pInt
|
||||
curl_fourier = (+tensorField_fourier(l,3,1,j,k)*xi(2,1,j,k)&
|
||||
|
@ -751,7 +751,7 @@ real(pReal) function utilities_curlRMS()
|
|||
-tensorField_fourier(l,1,1,j,k)*xi(2,1,j,k))*TWOPIIMG
|
||||
enddo
|
||||
utilities_curlRMS = utilities_curlRMS + &
|
||||
2.0_pReal*sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)
|
||||
sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)! this layer (DC) does not have a conjugate complex counterpart (if grid(1) /= 1)
|
||||
do l = 1_pInt, 3_pInt
|
||||
curl_fourier = (+tensorField_fourier(l,3,grid1Red,j,k)*xi(2,grid1Red,j,k)&
|
||||
-tensorField_fourier(l,2,grid1Red,j,k)*xi(3,grid1Red,j,k))*TWOPIIMG
|
||||
|
@ -761,8 +761,9 @@ real(pReal) function utilities_curlRMS()
|
|||
-tensorField_fourier(l,1,grid1Red,j,k)*xi(2,grid1Red,j,k))*TWOPIIMG
|
||||
enddo
|
||||
utilities_curlRMS = utilities_curlRMS + &
|
||||
2.0_pReal*sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)
|
||||
sum(real(curl_fourier)**2.0_pReal + aimag(curl_fourier)**2.0_pReal)! this layer (Nyquist) does not have a conjugate complex counterpart (if grid(1) /= 1)
|
||||
enddo; enddo
|
||||
|
||||
call MPI_Allreduce(MPI_IN_PLACE,utilities_curlRMS,1,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)
|
||||
utilities_curlRMS = sqrt(utilities_curlRMS) * wgt
|
||||
if(grid(1) == 1_pInt) utilities_curlRMS = utilities_curlRMS * 0.5_pReal ! counted twice in case of grid(1) == 1
|
||||
|
|
Loading…
Reference in New Issue