-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathPSEDirectSum.f90
596 lines (517 loc) · 22.1 KB
/
PSEDirectSum.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
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
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
module PSEDirectSumModule
use NumberKindsModule
use LoggerModule
use ParticlesModule
use EdgesModule, only : MaxEdgeLength
use PolyMesh2dModule
use FieldModule
use MPISetupModule
use SphereGeomModule, only : ChordDistance, SphereDistance, SphereProjection
implicit none
include 'mpif.h'
private
public PSE, New, Delete
public PSEPlaneInterpolateScalar
public PSEPlaneGradientAtParticles, PSEPlaneSecondPartialsAtParticles, PSEPlaneLaplacianAtParticles
public PSESphereInterpolateScalar, PSESphereGradientAtParticles, PSESphereDivergenceAtParticles, PSESphereLaplacianAtParticles
public bivariateFirstDerivativeKernel8, bivariateLaplacianKernel8
public PSEPlaneDoubleDotProductAtParticles
type PSE
real(kreal) :: eps
contains
final :: deletePrivate
end type
interface New
module procedure newPrivate
end interface
interface Delete
module procedure deletePrivate
end interface
!
!----------------
! Logging
!----------------
!
logical(klog), save :: logInit = .FALSE.
type(Logger) :: log
character(len=28), save :: logKey = 'PSE'
integer(kint), parameter :: logLevel = DEBUG_LOGGING_LEVEL
contains
!----------------
!
! Public functions
!
!----------------
subroutine newPrivate( self, aMesh, radiusMultiplier )
type(PSE), intent(out) :: self
type(PolyMesh2d), intent(in) :: aMesh
real(kreal), intent(in), optional :: radiusMultiplier
!
real(kreal) :: multiple
integer(kint) :: nP
if ( .NOT. logInit ) call InitLogger(log, procRank)
if ( present(radiusMultiplier) ) then
multiple = radiusMultiplier
else
multiple = 2.0_kreal
endif
self%eps = multiple * MaxEdgeLength( aMesh%edges, aMesh%particles)
end subroutine
subroutine deletePrivate(self)
type(PSE), intent(inout) :: self
end subroutine
pure function PSEPlaneInterpolateScalar(self, aMesh, scalarField, interpLoc )
real(kreal) :: PSEPlaneInterpolateScalar
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
real(kreal), intent(in) :: interpLoc(3)
!
integer(kint) :: j
real(kreal) :: kIn, xj(3)
PSEPlaneInterpolateScalar = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles,j)
kIn = ChordDistance(interpLoc,xj)/self%eps
PSEPlaneInterpolateScalar = PSEPlaneInterpolateScalar + scalarField%scalar(j) * aMesh%particles%area(j) * &
bivariateDeltaKernel8( kIn ) / (self%eps**2)
endif
enddo
end function
subroutine PSEPlaneGradientAtParticles( self, aMesh, scalarField, scalarGrad, particlesMPI )
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
type(Field), intent(inout) :: scalarGrad
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), kOut
call SetFieldToZero(scalarGrad)
scalarGrad%N = aMesh%particles%N
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logKey), " entering PSEGradientAtParticles...")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = ChordDistance(xi,xj)/self%eps
scalarGrad%xComp(i) = scalarGrad%xComp(i) + ( scalarField%scalar(j) + scalarField%scalar(i)) *&
(xi(1)-xj(1)) * bivariateFirstDerivativeKernel8(kIn) / (self%eps**3) * aMesh%particles%area(j)
scalarGrad%yComp(i) = scalarGrad%yComp(i) + ( scalarField%scalar(j) + scalarField%scalar(i)) *&
(xi(2)-xj(2)) * bivariateFirstDerivativeKernel8(kIn) / (self%eps**3) * aMesh%particles%area(j)
endif
enddo
enddo
do i = 0, numProcs - 1
call MPI_BCAST(scalarGrad%xComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
call MPI_BCAST(scalarGrad%yComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
enddo
call MultiplyFieldByScalar(scalarGrad, 1.0_kreal/self%eps)
end subroutine
subroutine PSESphereGradientAtParticles( self, amesh, scalarField, scalarGrad, particlesMPI )
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
type(Field), intent(inout) :: scalarGrad
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), grad(3), proj(3,3), kOut
call SetFieldToZero(scalarGrad)
scalarGrad%N = aMesh%particles%N
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logkey), " entering PSESphereGradientAtParticles.")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
proj = SphereProjection(xi)
grad = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord( aMesh%particles, j)
kIn = SphereDistance(xj, xi) / self%eps
kOut = bivariateFirstDerivativeKernel8( kIn ) / self%eps**2
grad(1) = grad(1) + ( scalarField%scalar(j) + scalarField%scalar(i)) * ( xi(1) - xj(1) ) * &
kOut * aMesh%particles%area(j)
grad(2) = grad(2) + ( scalarField%scalar(j) + scalarField%scalar(i)) * ( xi(2) - xj(2) ) * &
kOut * aMesh%particles%area(j)
grad(3) = grad(3) + ( scalarField%scalar(j) + scalarField%scalar(i)) * ( xi(3) - xj(3) ) * &
kOut * aMesh%particles%area(j)
endif
enddo
grad = MATMUL( proj, grad)
scalarGrad%xComp(i) = grad(1)
scalarGrad%yComp(i) = grad(2)
scalarGrad%zComp(i) = grad(3)
enddo
do i = 0, numProcs - 1
call MPI_BCAST( scalarGrad%xComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
call MPI_BCAST( scalarGrad%yComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
call MPI_BCAST( scalarGrad%zComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
enddo
call MultiplyFieldByScalar(scalarGrad, 1.0_kreal/(self%eps**2))
end subroutine
subroutine PSEPlaneSecondPartialsAtParticles( self, aMesh, scalarGrad, secondPartials, particlesMPI )
type(PSE), intent(in) :: self
class(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarGrad
type(Field), intent(inout) :: secondPartials
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), dxx, dxy, dyx, dyy, kOut
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logKey), " entering PSESecondPartialsAtParticles...")
call SetFieldToZero(secondPartials)
secondPartials%N = aMesh%particles%N
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord( aMesh%particles, i)
dxx = 0.0_kreal
dxy = 0.0_kreal
dyx = 0.0_kreal
dyy = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = ChordDistance(xi,xj)/self%eps
kOUt = bivariateFirstDerivativeKernel8(kIn) / self%eps**2
dxx = dxx + (scalarGrad%xComp(j) + scalarGrad%xComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
dxy = dxy + (scalarGrad%xComp(j) + scalarGrad%xComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
dyx = dyx + (scalarGrad%yComp(j) + scalarGrad%yComp(i)) * ( xi(1) - xj(1)) /self%eps * &
kOut * aMesh%particles%area(j)
dyy = dyy + (scalarGrad%yComp(j) + scalarGrad%yComp(i)) * ( xi(2) - xj(2)) /self%eps * &
kOut * aMesh%particles%area(j)
endif
enddo
secondPartials%xComp(i) = dxx
secondPartials%yComp(i) = 0.5_kreal * ( dxy + dyx )
secondPartials%zComp(i) = dyy
enddo
do i = 0, numProcs - 1
call MPI_BCAST(secondPartials%xComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCOde)
call MPI_BCAST(secondPartials%yComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCOde)
call MPI_BCAST(secondPartials%zComp(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCOde)
enddo
call MultiplyFieldByScalar( secondPartials, 1.0_kreal / self%eps )
end subroutine
subroutine PSEPlaneDoubleDotProductAtParticles( self, aMesh, vectorField, doubleDot, particlesMPI)
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: vectorField
type(Field), intent(inout) :: doubleDot
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), ux, uy, vx, vy, kOut
call SetFieldToZero(doubleDot)
doubleDot%N = aMesh%particles%N
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logKey), " entering PSEDoubleDotProduct ...")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
ux = 0.0_kreal
uy = 0.0_kreal
vx = 0.0_kreal
vy = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = ChordDistance(xi,xj)/self%eps
kOut = bivariateFirstDerivativeKernel8( kIn ) / self%eps**2
ux = ux + (vectorField%xComp(j) + vectorField%xComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
uy = uy + (vectorField%xComp(j) + vectorField%xComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
vx = vx + (vectorField%yComp(j) + vectorField%yComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
vy = vy + (vectorField%yComp(j) + vectorField%yComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
endif
enddo
doubleDot%scalar(i) = (ux * ux + 2.0_kreal * uy * vx + vy * vy)/self%eps**2
enddo
do i = 0, numProcs - 1
call MPI_BCAST( doubleDot%scalar(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode )
enddo
end subroutine
subroutine PSESphereDoubleDotProductAtParticles( self, aMesh, vectorField, doubleDot, particlesMPI)
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: vectorField
type(Field), intent(inout) :: doubleDot
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), ux, uy, uz, vx, vy, vz, wx, wy, wz, kOut
call SetFieldToZero(doubleDot)
doubleDot%N = aMesh%particles%N
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
ux = 0.0_kreal
uy = 0.0_kreal
uz = 0.0_kreal
vx = 0.0_kreal
vy = 0.0_kreal
vz = 0.0_kreal
wx = 0.0_kreal
wy = 0.0_kreal
wz = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = SphereDistance(xi,xj) / self%eps
kOut = bivariateFirstDerivativeKernel8( kIn ) / self%eps**2
ux = ux + (vectorField%xComp(j) + vectorField%xComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
uy = uy + (vectorField%xComp(j) + vectorField%xComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
uz = uz + (vectorField%xComp(j) + vectorField%xComp(i)) * ( xi(3) - xj(3) )/self%eps * &
kOut * aMesh%particles%area(j)
vx = vx + (vectorField%yComp(j) + vectorField%yComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
vy = vy + (vectorField%yComp(j) + vectorField%yComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
vz = vz + (vectorField%yComp(j) + vectorField%yComp(i)) * ( xi(3) - xj(3) )/self%eps * &
kOut * aMesh%particles%area(j)
wx = wx + (vectorField%zComp(j) + vectorField%yComp(i)) * ( xi(1) - xj(1) )/self%eps * &
kOut * aMesh%particles%area(j)
wy = wy + (vectorField%zComp(j) + vectorField%yComp(i)) * ( xi(2) - xj(2) )/self%eps * &
kOut * aMesh%particles%area(j)
wz = wz + (vectorField%zComp(j) + vectorField%yComp(i)) * ( xi(3) - xj(3) )/self%eps * &
kOut * aMesh%particles%area(j)
endif
enddo
doubleDot%scalar(i) = (ux*ux + vy*vy + wz*wz + 2.0_kreal*( uy*vx + uz*wx + vz*wy))/self%eps**2
enddo
end subroutine
subroutine PSEDoubleDotProductArrays( self, doubleDot, xIn, yIn, uIn, vIn, areaIn, activeMask, mpiParticles )
type(PSE), intent(in) :: self
real(kreal), intent(out) :: doubleDot(:)
real(kreal), intent(in) :: xIn(:)
real(kreal), intent(in) :: yIn(:)
real(kreal), intent(in) :: uIn(:)
real(kreal), intent(in) :: vIn(:)
real(kreal), intent(in) :: areaIn(:)
logical(klog), intent(in) :: activeMask(:)
type(MPISetup), intent(in) :: mpiParticles
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, ux, uy, vx, vy, kernel
doubleDot(mpiParticles%indexStart(procRank):mpiParticles%indexEnd(procRank)) = 0.0_kreal
do i = mpiParticles%indexStart(procRank), mpiParticles%indexEnd(procRank)
ux = 0.0_kreal
uy = 0.0_kreal
vx = 0.0_kreal
vy = 0.0_kreal
if ( activeMask(j) ) then
kIn = sqrt( (xIn(i) - xIn(j))**2 + (yIn(i) - yIn(j))**2 ) / self%eps
kernel = bivariateFirstDerivativeKernel8( kIn ) / self%eps**2
ux = ux + ( uIn(j) + uIn(i) ) * ( xIn(i) - xIn(j) ) / self%eps * kernel * areaIn(j)
uy = uy + ( uIn(j) + uIn(i) ) * ( yIn(i) - yIn(j) ) / self%eps * kernel * areaIn(j)
vx = vx + ( vIn(j) + vIn(i) ) * ( xIn(i) - xIn(j) ) / self%eps * kernel * areaIn(j)
vy = vy + ( vIn(j) + vIn(i) ) * ( yIn(i) - yIn(j) ) / self%eps * kernel * areaIn(j)
endif
doubleDot(i) = ( ux * ux + 2.0_kreal * uy * vx + vy * vy ) / self%eps**2
enddo
do i = 0, numProcs - 1
call MPI_BCAST( doubleDot(mpiParticles%indexStart(i):mpiParticles%indexEnd(i)), mpiParticles%messageLength(i), &
MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
enddo
end subroutine
subroutine PSEPlaneLaplacianAtParticles(self, aMesh, scalarField, scalarLap, particlesMPI )
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
type(Field), intent(inout) :: scalarLap
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3)
call SetFieldToZero(scalarLap)
scalarLap%N = aMesh%particles%N
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logKey), " entering PSELaplacianAtParticles...")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = ChordDistance(xi,xj)/self%eps
scalarLap%scalar(i) = scalarLap%scalar(i) + (scalarField%scalar(j) - scalarField%scalar(i)) * &
bivariateLaplacianKernel8(kIN)/(self%eps**2) * aMesh%particles%area(j)
endif
enddo
enddo
do i = 0, numProcs-1
call MPI_BCAST(scalarLap%scalar(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
enddo
call MultiplyFieldByScalar(scalarLap, 1.0_kreal/(self%eps**2))
end subroutine
subroutine PSESphereLaplacianAtParticles( self, aMesh, scalarField, scalarLap, particlesMPI )
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
type(Field), intent(inout) :: scalarLap
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3)
call SetFieldToZero(scalarLap)
scalarLap%N = aMesh%particles%N
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logKey), " entering PSELaplacianAtParticles...")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = SphereDistance(xi,xj)/self%eps
scalarLap%scalar(i) = scalarLap%scalar(i) + (scalarField%scalar(j) - scalarField%scalar(i)) * &
bivariateLaplacianKernel8(kIN)/(self%eps**2) * aMesh%particles%area(j)
endif
enddo
enddo
do i = 0, numProcs-1
call MPI_BCAST(scalarLap%scalar(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode)
enddo
call MultiplyFieldByScalar(scalarLap, 1.0_kreal/(self%eps**2))
end subroutine
pure function PSESphereInterpolateScalar( self, aMesh, scalarField, interpLoc)
real(kreal) :: PSESphereInterpolateScalar
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: scalarField
real(kreal), intent(in) :: interpLoc(3)
!
integer(kint) :: j
real(kreal) :: kIn, xj(3)
PSESphereInterpolateScalar = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = SphereDistance(xj, interpLoc)/self%eps
PSESphereInterpolateScalar = PSESphereInterpolateScalar + scalarField%scalar(j) * aMesh%particles%area(j) * &
bivariateDeltaKernel8( kIn ) / (self%eps**2)
endif
enddo
end function
subroutine PSESphereDivergenceAtParticles( self, aMesh, vectorField, divergence, particlesMPI )
type(PSE), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(Field), intent(in) :: vectorField
type(Field), intent(inout) :: divergence
type(MPISetup), intent(in) :: particlesMPI
!
integer(kint) :: i, j, mpiErrCode
real(kreal) :: kIn, xi(3), xj(3), grad(3), proj(3,3), denom
call SetFieldToZero(divergence)
divergence%N = aMesh%particles%N
denom = self%eps**3
call LogMessage(log, DEBUG_LOGGING_LEVEL, trim(logkey), " entering PSESphereDivergenceAtParticles.")
do i = particlesMPI%indexStart(procRank), particlesMPI%indexEnd(procRank)
xi = PhysCoord(aMesh%particles, i)
proj = SphereProjection(xi)
grad = 0.0_kreal
do j = 1, aMesh%particles%N
if ( aMesh%particles%isActive(j) ) then
xj = PhysCoord(aMesh%particles, j)
kIn = SphereDistance(xj, xi)/self%eps
grad(1) = (xi(1) - xj(1)) * bivariateFirstDerivativeKernel8( kIn ) / denom
grad(2) = (xi(2) - xj(2)) * bivariateFirstDerivativeKernel8( kIn ) / denom
grad(3) = (xi(3) - xj(3)) * bivariateFirstDerivativeKernel8( kIn ) / denom
grad = MATMUL( proj, grad )
divergence%scalar(i) = divergence%scalar(i) + ( grad(1) * (vectorField%xComp(j) + vectorField%xComp(i)) + &
grad(2) * (vectorField%yComp(j) + vectorField%yComp(i) ) +&
grad(3) * (vectorField%zComp(j) + vectorField%zComp(i) )) * aMesh%particles%area(j)
endif
enddo
enddo
do i = 0, numProcs - 1
call MPI_BCAST( divergence%scalar(particlesMPI%indexStart(i):particlesMPI%indexEnd(i)), &
particlesMPI%messageLength(i), MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, mpiErrCode )
enddo
call MultiplyFieldByScalar(divergence, 1.0_kreal/self%eps)
end subroutine
!----------------
!
! Private functions
!
!----------------
pure function integralKernelInput( xi, xj, eps, geomKind )
real(kreal) :: integralKernelInput
real(kreal), intent(in) :: xi(3), xj(3), eps
integer(kint), intent(in) :: geomKind
if ( geomKind == PLANAR_GEOM ) then
integralKernelInput = ChordDistance(xi,xj)/eps
elseif ( geomKind == SPHERE_GEOM ) then
integralKernelInput = SphereDistance(xi,xj)/eps
else
integralKernelInput = 0.0_kreal
endif
end function
pure function bivariateDeltaKernel8( radialDist )
real(kreal) :: bivariateDeltaKernel8
real(kreal), intent(in) :: radialDist
bivariateDeltaKernel8 = (4.0_kreal - 6.0_kreal * radialDist**2 + 2.0_kreal * radialDist**4 - &
radialDist**6 / 6.0_kreal) * exp(-radialDist * radialDist) / PI
end function
pure function bivariateLaplacianKernel8( radialDist )
real(kreal) :: bivariateLaplacianKernel8
real(kreal), intent(in) :: radialDist
bivariateLaplacianKernel8 = (40.0_kreal - 40.0_kreal * radialDist**2 + 10.0_kreal * radialDist**4 - &
2.0_kreal * radialDist**6/3.0_kreal) * exp( -radialDist*radialDist) / PI
end function
pure function bivariateFirstDerivativeKernel8( radialDist)
real(kreal) :: bivariateFirstDerivativeKernel8
real(kreal), intent(in) :: radialDist
bivariateFirstDerivativeKernel8 = (-20.0_kreal + 20.0_kreal * radialDist**2 - 5.0_kreal * radialDist**4 + &
radialDist**6/3.0_kreal) * exp( - radialDist*radialDist ) / PI
end function
pure function bivariateSecondDerivativeKernel8( radialDist )
real(kreal) :: bivariateSecondDerivativeKernel8
real(kreal), intent(in) :: radialDist
bivariateSecondDerivativeKernel8 = (24.0_kreal -16.0_kreal * radialDist**2 + 2.0_kreal * radialDist**4) * &
exp(- radialDist**2) / PI
end function
pure function bivariateMixedSecondDerivativeKernel8( radialDist, x, y )
real(kreal) :: bivariateMixedSecondDerivativeKernel8
real(kreal), intent(in) :: radialDist, x, y
bivariateMixedSecondDerivativeKernel8 = x * y * (160.0_kreal - 120.0_kreal * radialDist**2 + &
24.0_kreal * radialDist**4 - 4.0_kreal * radialDist**6/3.0_kreal )* exp(-radialDist**2)/PI
end function
pure function trivariateDeltaKernel8( radialDist )
real(kreal) :: trivariateDeltaKernel8
real(kreal), intent(in) :: radialDist
trivariateDeltaKernel8 = (6.5625_kreal - 7.875 * radialDist**2 + 2.25_kreal * radialDist**4 - &
radialDist**6/6.0_kreal) * exp( - radialDist**2 ) / 5.568327996831708_kreal
end function
pure function trivariateFirstDerivativeKernel8( radialDist, eps )
real(kreal) :: trivariateFirstDerivativeKernel8
real(kreal), intent(in) :: radialDist
real(kreal), intent(in) :: eps
trivariateFirstDerivativeKernel8 = (-7.875_kreal/eps**2 + 9.0_kreal * radialDist**2/eps**4 - radialDist**4/eps**4 - &
13.125_kreal - 15.75_kreal * radialDist**2 + 4.5_kreal * radialDist**4 - &
radialDist**6/3.0_kreal ) * exp(-radialDist**2)
end function
subroutine InitLogger(aLog,rank)
! Initialize a logger for this module and processor
type(Logger), intent(out) :: aLog
integer(kint), intent(in) :: rank
write(logKey,'(A,A,I0.2,A)') trim(logKey),'_',rank,' : '
if ( rank == 0 ) then
call New(aLog,logLevel)
else
call New(aLog,ERROR_LOGGING_LEVEL)
endif
logInit = .TRUE.
end subroutine
end module