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
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
|
! $Header$
#include "cctk.h"
module NoExcision_mod
implicit none
CCTK_REAL, parameter :: zero = 0.0
integer, parameter :: wp = kind(zero)
CCTK_REAL, parameter :: one = 1.0_wp
CCTK_REAL, parameter :: two = 2.0_wp
CCTK_INT :: infnorm_handle, sum_handle
CCTK_INT :: sym_selector
CCTK_INT :: loop_counter
CCTK_REAL, dimension(16) :: infnormresid
CCTK_REAL, dimension(16) :: delta_new, delta_old, alpha, beta, absres
CCTK_REAL, dimension(16) :: sumred, lsumred, infred, linfred
character(len=18), dimension(16) :: red_names = (/ 'NoExcision::redgxx', &
'NoExcision::redgxy', &
'NoExcision::redgxz', &
'NoExcision::redgyy', &
'NoExcision::redgyz', &
'NoExcision::redgzz', &
'NoExcision::redkxx', &
'NoExcision::redkxy', &
'NoExcision::redkxz', &
'NoExcision::redkyy', &
'NoExcision::redkyz', &
'NoExcision::redkzz', &
'NoExcision::red ', &
'NoExcision::redx ', &
'NoExcision::redy ', &
'NoExcision::redz ' /)
character(len=5), dimension(16) :: var_names = (/ ' gxx', ' gxy', &
' gxz', ' gyy', &
' gyz', ' gzz', &
' kxx', ' kxy', &
' kxz', ' kyy', &
' kyz', ' kzz', &
' alp', 'betax', &
'betay', 'betaz' /)
CCTK_INT :: nx, ny, nz
integer :: ii
logical, dimension(16) :: cont = (/ (.true., ii=1, 16) /)
contains
! This is the basic residual calculation routine...
subroutine residual (val, mask, res, si, order )
CCTK_REAL, dimension(:,:,:), intent(in) :: val
CCTK_INT, dimension(:,:,:), intent(in) :: mask
CCTK_REAL, dimension(:,:,:), intent(out) :: res
CCTK_REAL, intent(in) :: si
CCTK_INT, intent(in) :: order
CCTK_REAL, parameter :: six = 6.0_wp
CCTK_REAL, parameter :: sixth = 1.0_wp / six
CCTK_REAL, parameter :: four = 4.0_wp
CCTK_REAL, parameter :: eighteen = 18.0_wp
CCTK_REAL, parameter :: eighteenth = 1.0_wp / eighteen
CCTK_REAL, parameter :: fifteen = 15.0_wp
CCTK_REAL, parameter :: sixty = 60.0_wp
CCTK_REAL, parameter :: sixtieth = 1.0_wp / sixty
integer :: i, j, k
select case (order)
case (2)
!$OMP PARALLEL DO PRIVATE(k,j,i) SCHEDULE(GUIDED)
do k = 2, nz-1
do j = 2, ny-1
do i = 2, nx-1
if (mask(i,j,k) > 0) then
res(i,j,k) = val(i+1,j,k) + val(i-1,j,k) + val(i,j+1,k) + &
val(i,j-1,k) + val(i,j,k+1) + val(i,j,k-1) - &
six * val(i,j,k)
res(i,j,k) = sign(one,si) * res(i,j,k)
end if
end do
end do
end do
!$OMP END PARALLEL DO
case (4)
!$OMP PARALLEL DO PRIVATE(k,j,i) SCHEDULE(GUIDED)
do k = 3, nz-2
do j = 3, ny-2
do i = 3, nx-2
if (mask(i,j,k) > 0) then
res(i,j,k) = - ( val(i+2,j,k) + val(i-2,j,k) + &
val(i,j+2,k) + val(i,j-2,k) + &
val(i,j,k+2) + val(i,j,k-2) ) + &
four * ( val(i+1,j,k) + val(i-1,j,k) + &
val(i,j+1,k) + val(i,j-1,k) + &
val(i,j,k+1) + val(i,j,k-1) ) - &
eighteen * val(i,j,k)
res(i,j,k) = sign(one,si) * res(i,j,k)
end if
end do
end do
end do
!$OMP END PARALLEL DO
case (6)
!$OMP PARALLEL DO PRIVATE(k,j,i) SCHEDULE(GUIDED)
do k = 4, nz-3
do j = 4, ny-3
do i = 4, nx-3
if (mask(i,j,k) > 0) then
res(i,j,k) = ( val(i+3,j,k) + val(i-3,j,k) + &
val(i,j+3,k) + val(i,j-3,k) + &
val(i,j,k+3) + val(i,j,k-3) ) - &
six * ( val(i+2,j,k) + val(i-2,j,k) + &
val(i,j+2,k) + val(i,j-2,k) + &
val(i,j,k+2) + val(i,j,k-2) ) + &
fifteen * ( val(i+1,j,k) + val(i-1,j,k) + &
val(i,j+1,k) + val(i,j-1,k) + &
val(i,j,k+1) + val(i,j,k-1) ) - &
sixty * val(i,j,k)
res(i,j,k) = sign(one,si) * res(i,j,k)
end if
end do
end do
end do
!$OMP END PARALLEL DO
case default
call CCTK_WARN ( 0, 'Internal error: Order out of range' )
end select
end subroutine residual
! This is the routine calculating residuals for all variables, that
! have not yet converged.
subroutine residual_all ( v1, v2, v3, v4, v5, v6, v7, v8, &
v9, v10, v11, v12, v13, v14, v15, v16, &
r1, r2, r3, r4, r5, r6, r7, r8, &
r9, r10, r11, r12, r13, r14, r15, r16, &
mask, si, order )
CCTK_REAL, dimension(:,:,:), intent(in) :: v1, v2, v3, v4, v5, v6, &
v7, v8, v9, v10, v11, v12, &
v13, v14, v15, v16
CCTK_REAL, dimension(:,:,:), intent(out) :: r1, r2, r3, r4, r5, r6, &
r7, r8, r9, r10, r11, r12, &
r13, r14, r15, r16
CCTK_INT, dimension(:,:,:), intent(in) :: mask
CCTK_REAL, intent(in) :: si
CCTK_INT, intent(in) :: order
if ( cont(1) ) call residual ( v1, mask, r1, si, order )
if ( cont(2) ) call residual ( v2, mask, r2, si, order )
if ( cont(3) ) call residual ( v3, mask, r3, si, order )
if ( cont(4) ) call residual ( v4, mask, r4, si, order )
if ( cont(5) ) call residual ( v5, mask, r5, si, order )
if ( cont(6) ) call residual ( v6, mask, r6, si, order )
if ( cont(7) ) call residual ( v7, mask, r7, si, order )
if ( cont(8) ) call residual ( v8, mask, r8, si, order )
if ( cont(9) ) call residual ( v9, mask, r9, si, order )
if ( cont(10) ) call residual ( v10, mask, r10, si, order )
if ( cont(11) ) call residual ( v11, mask, r11, si, order )
if ( cont(12) ) call residual ( v12, mask, r12, si, order )
if ( cont(13) ) call residual ( v13, mask, r13, si, order )
if ( cont(14) ) call residual ( v14, mask, r14, si, order )
if ( cont(15) ) call residual ( v15, mask, r15, si, order )
if ( cont(16) ) call residual ( v16, mask, r16, si, order )
end subroutine residual_all
! This routine multiplies two variables. Only does it according
! to the mask and only for variables that have not converged yet.
subroutine multiply ( u1, u2, u3, u4, u5, u6, u7, u8, &
u9, u10, u11, u12, u13, u14, u15, u16, &
v1, v2, v3, v4, v5, v6, v7, v8, &
v9, v10, v11, v12, v13, v14, v15, v16, &
r1, r2, r3, r4, r5, r6, r7, r8, &
r9, r10, r11, r12, r13, r14, r15, r16, &
mask, weight, do_inf_reduction, order )
CCTK_REAL, dimension(:,:,:), intent(in) :: u1, u2, u3, u4, u5, u6, &
u7, u8, u9, u10, u11, u12, &
u13, u14, u15, u16
CCTK_REAL, dimension(:,:,:), intent(in) :: v1, v2, v3, v4, v5, v6, &
v7, v8, v9, v10, v11, v12, &
v13, v14, v15, v16
CCTK_REAL, dimension(:,:,:), intent(out) :: r1, r2, r3, r4, r5, r6, &
r7, r8, r9, r10, r11, r12, &
r13, r14, r15, r16
CCTK_INT, dimension(:,:,:), intent(in) :: mask
CCTK_REAL, dimension(:,:,:), intent(in) :: weight
LOGICAL, intent(in) :: do_inf_reduction
CCTK_INT, intent(in) :: order
CCTK_INT i, j, k, offset
offset = order / 2
lsumred = 0.0
if (do_inf_reduction) linfred = 0.0
!$OMP PARALLEL DO SCHEDULE(guided) PRIVATE(k,j,i) REDUCTION(+:lsumred) &
!$OMP REDUCTION(max:linfred)
do k = 1 + offset, nz - offset
do j = 1 + offset, ny - offset
do i = 1 + offset, nx - offset
if ( mask(i,j,k) > 0 ) then
if ( cont(1) ) then
r1(i,j,k) = u1(i,j,k)*v1(i,j,k)
lsumred(1) = lsumred(1) + weight(i,j,k)*r1(i,j,k)
if (do_inf_reduction) linfred(1) = max(linfred(1),r1(i,j,k))
end if
if ( cont(2) ) then
r2(i,j,k) = u2(i,j,k)*v2(i,j,k)
lsumred(2) = lsumred(2) + weight(i,j,k)*r2(i,j,k)
if (do_inf_reduction) linfred(2) = max(linfred(2),r2(i,j,k))
end if
if ( cont(3) ) then
r3(i,j,k) = u3(i,j,k)*v3(i,j,k)
lsumred(3) = lsumred(3) + weight(i,j,k)*r3(i,j,k)
if (do_inf_reduction) linfred(3) = max(linfred(3),r3(i,j,k))
end if
if ( cont(4) ) then
r4(i,j,k) = u4(i,j,k)*v4(i,j,k)
lsumred(4) = lsumred(4) + weight(i,j,k)*r4(i,j,k)
if (do_inf_reduction) linfred(4) = max(linfred(4),r4(i,j,k))
end if
if ( cont(5) ) then
r5(i,j,k) = u5(i,j,k)*v5(i,j,k)
lsumred(5) = lsumred(5) + weight(i,j,k)*r5(i,j,k)
if (do_inf_reduction) linfred(5) = max(linfred(5),r5(i,j,k))
end if
if ( cont(6) ) then
r6(i,j,k) = u6(i,j,k)*v6(i,j,k)
lsumred(6) = lsumred(6) + weight(i,j,k)*r6(i,j,k)
if (do_inf_reduction) linfred(6) = max(linfred(6),r6(i,j,k))
end if
if ( cont(7) ) then
r7(i,j,k) = u7(i,j,k)*v7(i,j,k)
lsumred(7) = lsumred(7) + weight(i,j,k)*r7(i,j,k)
if (do_inf_reduction) linfred(7) = max(linfred(7),r7(i,j,k))
end if
if ( cont(8) ) then
r8(i,j,k) = u8(i,j,k)*v8(i,j,k)
lsumred(8) = lsumred(8) + weight(i,j,k)*r8(i,j,k)
if (do_inf_reduction) linfred(8) = max(linfred(8),r8(i,j,k))
end if
if ( cont(9) ) then
r9(i,j,k) = u9(i,j,k)*v9(i,j,k)
lsumred(9) = lsumred(9) + weight(i,j,k)*r9(i,j,k)
if (do_inf_reduction) linfred(9) = max(linfred(9),r9(i,j,k))
end if
if ( cont(10) ) then
r10(i,j,k) = u10(i,j,k)*v10(i,j,k)
lsumred(10) = lsumred(10) + weight(i,j,k)*r10(i,j,k)
if (do_inf_reduction) linfred(10) = max(linfred(10),r10(i,j,k))
end if
if ( cont(11) ) then
r11(i,j,k) = u11(i,j,k)*v11(i,j,k)
lsumred(11) = lsumred(11) + weight(i,j,k)*r11(i,j,k)
if (do_inf_reduction) linfred(11) = max(linfred(11),r11(i,j,k))
end if
if ( cont(12) ) then
r12(i,j,k) = u12(i,j,k)*v12(i,j,k)
lsumred(12) = lsumred(12) + weight(i,j,k)*r12(i,j,k)
if (do_inf_reduction) linfred(12) = max(linfred(12),r12(i,j,k))
end if
if ( cont(13) ) then
r13(i,j,k) = u13(i,j,k)*v13(i,j,k)
lsumred(13) = lsumred(13) + weight(i,j,k)*r13(i,j,k)
if (do_inf_reduction) linfred(13) = max(linfred(13),r13(i,j,k))
end if
if ( cont(14) ) then
r14(i,j,k) = u14(i,j,k)*v14(i,j,k)
lsumred(14) = lsumred(14) + weight(i,j,k)*r14(i,j,k)
if (do_inf_reduction) linfred(14) = max(linfred(14),r14(i,j,k))
end if
if ( cont(15) ) then
r15(i,j,k) = u15(i,j,k)*v15(i,j,k)
lsumred(15) = lsumred(15) + weight(i,j,k)*r15(i,j,k)
if (do_inf_reduction) linfred(15) = max(linfred(15),r15(i,j,k))
end if
if ( cont(16) ) then
r16(i,j,k) = u16(i,j,k)*v16(i,j,k)
lsumred(16) = lsumred(16) + weight(i,j,k)*r16(i,j,k)
if (do_inf_reduction) linfred(16) = max(linfred(16),r16(i,j,k))
end if
end if
end do
end do
end do
!$OMP END PARALLEL DO
end subroutine multiply
! This routine takes the product of a scalar and a variable and
! adds it to the product of another scalar and variable. Does it
! according to the mask and only for variables that have not yet
! converged.
subroutine multiply_sum ( u1, u2, u3, u4, u5, u6, u7, u8, &
u9, u10, u11, u12, u13, u14, u15, u16, &
v1, v2, v3, v4, v5, v6, v7, v8, &
v9, v10, v11, v12, v13, v14, v15, v16, &
c1, c2, mask )
CCTK_REAL, dimension(:,:,:), intent(inout) :: u1, u2, u3, u4, u5, u6, &
u7, u8, u9, u10, u11, u12, &
u13, u14, u15, u16
CCTK_REAL, dimension(:,:,:), intent(in) :: v1, v2, v3, v4, v5, v6, &
v7, v8, v9, v10, v11, v12, &
v13, v14, v15, v16
CCTK_REAL, dimension(16), intent(in) :: c1, c2
CCTK_INT, dimension(:,:,:), intent(in) :: mask
!$OMP PARALLEL
if ( cont(1) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u1 = c1(1)*u1 + c2(1)*v1
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(2) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u2 = c1(2)*u2 + c2(2)*v2
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(3) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u3 = c1(3)*u3 + c2(3)*v3
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(4) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u4 = c1(4)*u4 + c2(4)*v4
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(5) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u5 = c1(5)*u5 + c2(5)*v5
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(6) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u6 = c1(6)*u6 + c2(6)*v6
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(7) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u7 = c1(7)*u7 + c2(7)*v7
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(8) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u8 = c1(8)*u8 + c2(8)*v8
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(9) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u9 = c1(9)*u9 + c2(9)*v9
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(10) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u10 = c1(10)*u10 + c2(10)*v10
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(11) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u11 = c1(11)*u11 + c2(11)*v11
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(12) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u12 = c1(12)*u12 + c2(12)*v12
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(13) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u13 = c1(13)*u13 + c2(13)*v13
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(14) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u14 = c1(14)*u14 + c2(14)*v14
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(15) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u15 = c1(15)*u15 + c2(15)*v15
!$OMP END WORKSHARE NOWAIT
endif
if ( cont(16) ) then
!$OMP WORKSHARE
where ( mask > 0 ) u16 = c1(16)*u16 + c2(16)*v16
!$OMP END WORKSHARE NOWAIT
endif
!$OMP END PARALLEL
end subroutine multiply_sum
end module NoExcision_mod
|