-
Notifications
You must be signed in to change notification settings - Fork 12
/
Copy pathm_force.f90
185 lines (135 loc) · 4.95 KB
/
m_force.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
module m_force
use m_parameters, only : kfmax, famp, force_type
integer*4 :: n_forced_nodes, n_forced_nodes_total
! coordinated of the forced nodes
integer, allocatable :: ifn(:), jfn(:), kfn(:), k_shell(:)
contains
subroutine m_force_init
use m_parameters
use x_fftw
implicit none
integer :: i, j, k, n, n_shell
! if flow is not forced, return
if (flow_type .ne. 1) return
if (task.ne.'hydro') return
select case (force_type)
case (1)
! Machiels forcing (see article in PRL #79(18) p.3411)
write(out,*) "Forcing #1: Machiels forcing - setting up"
call flush(out)
! find out how many nodes are we forcing and book them
n_forced_nodes = 0
n_forced_nodes_total = 0
do k = 1,nz
do j = 1,ny
do i = 1,nx
n_shell = nint(sqrt(real(akx(i)**2 + aky(k)**2 + akz(j)**2, 4)))
if (n_shell .gt. 0 .and. n_shell .le. kfmax) then
n_forced_nodes = n_forced_nodes + 1
end if
end do
end do
end do
! reducing to the master process to find out the total number of forced nodes
!!$ write(out,*) 'before reducing'; call flush(out)
count = 1
call MPI_REDUCE(n_forced_nodes, n_forced_nodes_total, count, &
MPI_INTEGER4,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
! writing out the # of forced nodes
write(out,*) 'Number of forced nodes for this process:',n_forced_nodes
if (myid.eq.0) write(out,*) ' total number:',n_forced_nodes_total
call flush(out)
! allocating arrays for the coordinates of the forced nodes
allocate(ifn(n_forced_nodes), jfn(n_forced_nodes), kfn(n_forced_nodes), &
k_shell(n_forced_nodes))
! filling up the arrays
n = 0
do k = 1,nz
do j = 1,ny
do i = 1,nx
n_shell = nint(sqrt(real(akx(i)**2 + aky(k)**2 + akz(j)**2, 4)))
if (n_shell .gt. 0 .and. n_shell .le. kfmax) then
n = n + 1
ifn(n) = i
jfn(n) = j
kfn(n) = k
k_shell(n) = n_shell
end if
end do
end do
end do
!!$ ! writing out the nodes
!!$ do n = 1,n_forced_nodes
!!$ write(out,"(3i4)") ifn(n),jfn(n),kfn(n)
!!$ end do
!!$ call flush(out)
case default
write(out,*) 'WRONG FORCE TYPE:',force_type
write(out,*) 'STOPPING'
call flush(out)
stop
end select
return
end subroutine m_force_init
!================================================================================
subroutine force_velocity
! adding forcing to the arrays wrk(:,:,:,1:3) that already contain the RHS for velocities
use m_openmpi
use m_io
use m_parameters
use m_fields
use m_work
use x_fftw
use m_stats
implicit none
integer :: i, j, k, n_shell, n
real*8 :: fac, fac2
select case (force_type)
case (1)
! Machiels forcing (see article in PRL #79(18) p.3411)
! write(out,*) "Machiels forcing"; call flush(out)
e_spec = zip
e_spec1 = zip
hits = 0
hits1 = 0
! need this normalization factor because the FFT is unnormalized
fac = one / real(nx*ny*nz_all)**2
! assembling the total energy in each shell
do k = 1,nz
do j = 1,ny
do i = 1,nx
n_shell = nint(sqrt(real(akx(i)**2 + aky(k)**2 + akz(j)**2, 4)))
if (n_shell .gt. 0 .and. n_shell .le. kfmax) then
fac2 = fac * (fields(i,j,k,1)**2 + fields(i,j,k,2)**2 + fields(i,j,k,3)**2)
if (akx(i).eq.0.d0) fac2 = fac2 * 0.5d0
e_spec1(n_shell) = e_spec1(n_shell) + fac2
end if
end do
end do
end do
! reducing the number of hits and energy to two arrays on master node
count = kfmax
call MPI_REDUCE(e_spec1,e_spec,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
! getting the total energy in the region [0:kfmax] by integrating the spectrum
if (myid.eq.0) energy = sum(e_spec(1:kfmax))
! broadcasting the current energy in the forcing range of wavenumbers
count = 1
call MPI_BCAST(energy,count,MPI_REAL8,0,MPI_COMM_TASK,mpi_err)
! now applying the forcing to the RHS for velocities (wrk(:,:,:,1:3))
fac = FAMP / energy
do n = 1,n_forced_nodes
n_shell = k_shell(n)
i = ifn(n)
j = jfn(n)
k = kfn(n)
wrk(i,j,k,1) = wrk(i,j,k,1) + fac * fields(i,j,k,1)
wrk(i,j,k,2) = wrk(i,j,k,2) + fac * fields(i,j,k,2)
wrk(i,j,k,3) = wrk(i,j,k,3) + fac * fields(i,j,k,3)
end do
case default
write(out,*) "WRONG FORCE_TYPE:",force_type
stop
end select
return
end subroutine force_velocity
end module m_force