-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmatrix.f95
245 lines (191 loc) · 6.1 KB
/
matrix.f95
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
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
module matrix
use grid
use mpi_interface
contains
subroutine CalcVelocity(phi,psi,u,v,M,N,L)
implicit none
real*8, intent(in), dimension(M,N+2,L+2) :: phi,psi
real*8, intent(out), dimension(M,N,L) :: u,v
integer, intent(in) :: M,N,L
real*8, dimension(M,N,L) :: phix,phiy,psix,psiy
call Grad(phi,phix,phiy,M,N,L)
call Grad(psi,psix,psiy,M,N,L)
u=-phix-psiy
v=-phiy+psix
end subroutine CalcVelocity
subroutine Grad(p,px,py,M,N,L)
real*8, dimension(M,N+2,L+2), intent(in) :: p
real*8, dimension(M,N,L), intent(out) :: px,py
integer, intent(in) :: M,N,L
integer :: i
do i=2,N+1
px(:,i-1,:)=(0.5/hx)*(p(:,i+1,2:L+1)-p(:,i-1,2:L+1))
enddo
do i=2,L+1
py(:,:,i-1)=(0.5/hy)*(p(:,2:N+1,i+1)-p(:,2:N+1,i-1))
enddo
end subroutine Grad
subroutine MatMultForcing(p,q,M,N,L,hx,hy,hz,opt)
implicit none
real*8, dimension(M+2,N+2,L+2), intent(in) :: q !this is what gets INPUT. Layer of ghost cells around q.
real*8, dimension(M,N,L), intent(out) :: p !this is what gets OUTPUT. No ghost cells around p.
real*8, intent(in) :: hx,hz,hy
integer, intent(in) :: M,N,L,opt
integer :: i,j,k
real*8 :: pxx,pyy,pzz,pyz,theta,alpha,beta,px,py,pz,p1y,p2y,p3y
alpha=dbeta/dalpha!d2/d1
beta=dalpha*dbeta
if(opt==1)then !w equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
pxx=(1.0/hx**2)*(q(i,j+1,k)-2.0*q(i,j,k)+q(i,j-1,k))
pyy=(1.0/hy**2)*(q(i,j,k+1)-2.0*q(i,j,k)+q(i,j,k-1))
pzz=(1.0/hz**2)*(q(i+1,j,k)-2.0*q(i,j,k)+q(i-1,j,k))
p(i-1,j-1,k-1)=(pxx+pyy)/(1.0+beta)
enddo
enddo
enddo
elseif(opt==2)then !b equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
pxx=(1.0/hx**2)*(q(i,j+1,k)-2.0*q(i,j,k)+q(i,j-1,k))
pyy=(1.0/hy**2)*(q(i,j,k+1)-2.0*q(i,j,k)+q(i,j,k-1))
if(i==2)then
pzz=(1.0/hz**2)*(2.0*q(i,j,k)-5.0*q(i+1,j,k)+4.0*q(i+2,j,k)-q(i+3,j,k))
p1y=(0.5/hy)*(q(i,j,k+1)-q(i,j,k-1))
p2y=(0.5/hy)*(q(i+1,j,k+1)-q(i+1,j,k-1))
p3y=(0.5/hy)*(q(i+2,j,k+1)-q(i+2,j,k-1))
pyz=(1.0/hz)*(-1.5*p1y+2.0*p2y-0.5*p3y)
elseif(i==M+1)then
pzz=(1.0/hz**2)*(-q(i-3,j,k)+4.0*q(i-2,j,k)-5.0*q(i-1,j,k)+2.0*q(i,j,k))
p1y=(0.5/hy)*(q(i,j,k+1)-q(i,j,k-1))
p2y=(0.5/hy)*(q(i-1,j,k+1)-q(i-1,j,k-1))
p3y=(0.5/hy)*(q(i-2,j,k+1)-q(i-2,j,k-1))
pyz=(1.0/hz)*(0.5*p3y-2.0*p2y+1.5*p1y)
else
pzz=(1.0/hz**2)*(q(i+1,j,k)-2.0*q(i,j,k)+q(i-1,j,k))
pyz=(1.0/(4.0*hy*hz))*(q(i+1,j,k+1)-q(i+1,j,k-1)-q(i-1,j,k+1)+q(i-1,j,k-1))
endif
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=(dalpha/(1.0+beta))*(pxx+pyy+pzz) &
+(1.0/(dalpha*(1.0+beta)))*((cos(theta)**2)*pyy+2.0*cos(theta)*sin(theta)*pyz+(sin(theta)**2)*pzz)
enddo
enddo
enddo
elseif(opt==3)then !p equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
px=(q(i,j+1,k)-q(i,j-1,k))/(2.0*hx)
py=(q(i,j,k+1)-q(i,j,k-1))/(2.0*hy)
pz=(q(i+1,j,k)-q(i-1,j,k))/(2.0*hz)
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=(1.0/(1.0+beta))*((((sin(theta)**2)/(dalpha))+dalpha)*pz + (1.0/dalpha)*sin(theta)*cos(theta)*py - cos(theta)*px)
enddo
enddo
enddo
elseif(opt==4)then !phi equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
px=(q(i,j+1,k)-q(i,j-1,k))/(2.0*hx)
py=(q(i,j,k+1)-q(i,j,k-1))/(2.0*hy)
pz=(q(i+1,j,k)-q(i-1,j,k))/(2.0*hz)
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=(1.0/(1.0+beta))*pz
enddo
enddo
enddo
elseif(opt==5)then !psi equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
px=(q(i,j+1,k)-q(i,j-1,k))/(2.0*hx)
py=(q(i,j,k+1)-q(i,j,k-1))/(2.0*hy)
pz=(q(i+1,j,k)-q(i-1,j,k))/(2.0*hz)
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=(1.0/(dalpha*(1.0+beta)))*(sin(theta)*pz+cos(theta)*py)
enddo
enddo
enddo
elseif(opt==6)then !u equation
do j=2,N+1
do i=2,M+1
do k=2,L+1
px=(q(i,j+1,k)-q(i,j-1,k))/(2.0*hx)
py=(q(i,j,k+1)-q(i,j,k-1))/(2.0*hy)
pz=(q(i+1,j,k)-q(i-1,j,k))/(2.0*hz)
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=(1.0/(dalpha*(1.0+beta)))*(sin(theta)*pz+cos(theta)*py)
enddo
enddo
enddo
endif
end subroutine MatMultForcing
subroutine MatMult(p,q,M,N,L,hx,hy,hz)
implicit none
real*8, dimension(M+2,N+2,L+2), intent(in) :: q !this is what gets INPUT. Layer of ghost cells around q.
real*8, dimension(M,N,L), intent(out) :: p !this is what gets OUTPUT. No ghost cells around p.
real*8, intent(in) :: hx,hz,hy
integer, intent(in) :: M,N,L
integer :: i,j,k
real*8 :: pxx,pyy,pzz,pyz,theta,alpha,beta
alpha=dbeta/dalpha
beta=dalpha*dbeta
do j=2,N+1
do i=2,M+1
do k=2,L+1
pxx=(1.0/hx**2)*(q(i,j+1,k)-2.0*q(i,j,k)+q(i,j-1,k))
pyy=(1.0/hy**2)*(q(i,j,k+1)-2.0*q(i,j,k)+q(i,j,k-1))
pzz=(1.0/hz**2)*(q(i+1,j,k)-2.0*q(i,j,k)+q(i-1,j,k))
pyz=(1.0/(4.0*hy*hz))*(q(i+1,j,k+1)-q(i+1,j,k-1)-q(i-1,j,k+1)+q(i-1,j,k-1))
theta=(4.0*atan(1.0)/180.0)*dtheta
p(i-1,j-1,k-1)=pxx+pyy &
+(dalpha/(1.0+beta))*((cos(theta)**2)*pyy+2.0*cos(theta)*sin(theta)*pyz+(sin(theta)**2)*pzz) &
+(beta/(1.0+beta))*pzz
enddo
enddo
enddo
end subroutine MatMult
subroutine ConjGrad(x,f,tol,N,M,L,hx,hy,hz)
implicit none
real*8, dimension(M+2,N+2,L+2), intent(out) :: x !This is what gets output
real*8, dimension(M,N,L), intent(in) :: f !This is the right hand side
real*8, intent(in) :: tol
integer :: i,j,k
real*8, dimension(M,N,L) :: r,q,Ax
integer, intent(in) :: N,M,L
real*8 :: hx,hz,hy
!x(2:M+1,2:N+1,2:L+1) = f
r=f
rhs_sum=sum(r)
call mpi_allreduce(rhs_sum,rhs_sum_tot,1,mpi_real8,mpi_sum,comm2d,ierr)
!r=r-(rhs_sum_tot/(N_x*N_vert*N_y))
p(2:M+1,2:N+1,2:L+1)=r
rold=sum(r*r)
call mpi_allreduce(rold,rold_tot,1,mpi_real8,mpi_sum,comm2d,ierr)
do k=1,N_x**2
call conj_exchange(p,k)
p(1,:,:)=-p(2,:,:)
p(M+2,:,:)=-p(M+1,:,:)
!write(*,*) size(p,1),size(p,2),size(p,3)
call matmult(q,p,M,N,L,hx,hy,hz)
gam=sum(q*p(2:M+1,2:N+1,2:L+1))
call mpi_allreduce(gam,gam_tot,1,mpi_real8,mpi_sum,comm2d,ierr)
alpha_tot=rold_tot/gam_tot
x(2:M+1,2:N+1,2:L+1) = x(2:M+1,2:N+1,2:L+1)+alpha_tot*p(2:M+1,2:N+1,2:L+1)
r=r-alpha_tot*q
rnew=sum(r*r)
call mpi_allreduce(rnew,rnew_tot,1,mpi_real8,mpi_sum,comm2d,ierr)
if(sqrt(rnew_tot)<tol)then
exit
endif
p(2:M+1,2:N+1,2:L+1) = r + (rnew_tot/rold_tot)*p(2:M+1,2:N+1,2:L+1)
rold_tot=rnew_tot
write(*,*) k,rnew_tot
enddo
write(*,*) k, rnew
end subroutine ConjGrad
end module matrix