-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathssrfpackInterface.f90
391 lines (333 loc) · 14.8 KB
/
ssrfpackInterface.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
module SSRFPACKInterfaceModule
use NumberKindsModule
use OutputWriterModule
use LoggerModule
use ParticlesModule
use EdgesModule
use FacesModule
use PolyMesh2dModule
use FieldModule
use MPISetupModule
use SphereGeomModule
implicit none
include 'mpif.h'
private
!
!----------------
! Module types and public declarations
!----------------
!
public SSRFPACKInterface, DelaunayTriangulation, New, Delete
public SetScalarSourceData, SetVectorSourceData, SetSourceLagrangianParameter
public InterpolateScalar, InterpolateVector, InterpolateLagParam
public InterpolateScalarToUnifLatLonGrid, InterpolateVectorToUnifLatLonGrid
type DelaunayTriangulation
integer(kint), dimension(:), pointer :: list => null()
integer(kint), dimension(:), pointer :: lptr => null()
integer(kint), dimension(:), pointer :: lend => null()
contains
final :: deleteDelTri
end type
type SSRFPACKInterface
real(kreal), dimension(:,:), pointer :: grad1 => null()
real(kreal), dimension(:,:), pointer :: grad2 => null()
real(kreal), dimension(:,:), pointer :: grad3 => null()
real(kreal), dimension(:), pointer :: sigma1 => null()
real(kreal), dimension(:), pointer :: sigma2 => null()
real(kreal), dimension(:), pointer :: sigma3 => null()
contains
final :: deletePrivate
end type
integer(kint), parameter :: SIGMA_FLAG = 1
integer(kint), parameter :: GRAD_FLAG = 1
integer(kint), save :: startTriangle = 1
real(kreal), save :: sigmaTol = 0.01_kreal
!
!----------------
! Module interfaces
!----------------
!
interface New
module procedure newPrivate
module procedure newDelTri
end interface
interface Delete
module procedure deletePrivate
module procedure deleteDelTri
end interface
!
!----------------
! Logging
!----------------
!
logical(klog), save :: logInit = .FALSE.
type(Logger) :: log
character(len=28), save :: logKey = 'SSRFPACK'
integer(kint), parameter :: logLevel = DEBUG_LOGGING_LEVEL
character(len=MAX_STRING_LENGTH) :: logString
contains
!
!----------------
! public methods
!----------------
!
subroutine newPrivate(self, aMesh, nDim )
type(SSRFPACKInterface), intent(out) :: self
type(PolyMesh2d), intent(inout) :: aMesh
integer(kint), intent(in) :: nDim
if ( .NOT. logInit) call InitLogger(log, procRank)
allocate(self%grad1(3, aMesh%particles%N))
allocate(self%sigma1( 6 * aMesh%particles%N - 12 ))
if ( nDim == 3 ) then
allocate(self%grad2(3, aMesh%particles%N))
allocate(self%grad3(3, aMesh%particles%N))
allocate(self%sigma2( 6 * aMesh%particles%N - 12 ))
allocate(self%sigma3( 6 * aMesh%particles%N - 12 ))
endif
end subroutine
subroutine newDelTri(self, aMesh)
type(DelaunayTriangulation), intent(out) :: self
type(PolyMesh2d), intent(inout) :: aMesh
if ( .NOT. logInit) call InitLogger(log, procRank)
call ProjectParticlesToSphere(aMesh, 1.0_kreal)
allocate(self%list( 6 * aMesh%particles%N - 12 ))
allocate(self%lptr( 6 * amesh%particles%N - 12 ))
allocate(self%lend( aMesh%particles%N ))
call BuildDelaunayTriangulation(self, aMesh)
end subroutine
subroutine deleteDelTri(self)
type(DelaunayTriangulation), intent(inout) :: self
if ( associated(self%list)) then
deallocate(self%list)
deallocate(self%lptr)
deallocate(self%lend)
endif
end subroutine
subroutine deletePrivate( self )
type(SSRFPACKInterface), intent(inout) :: self
if ( associated(self%grad1) ) then
deallocate(self%grad1)
deallocate(self%sigma1)
endif
if ( associated(self%grad2) ) then
deallocate(self%grad2)
deallocate(self%grad3)
deallocate(self%sigma2)
deallocate(self%sigma3)
endif
end subroutine
subroutine SetVectorSourceData(self, aMesh, delTri, vectorField)
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: vectorField
integer(kint) :: i, errCode
real(kreal) :: dSig
do i = 1, aMesh%particles%N
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%xComp, delTri%list, delTri%lptr, delTri%lend, self%grad1(:,i), errCode)
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%yComp, delTri%list, delTri%lptr, delTri%lend, self%grad2(:,i), errCode)
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%zComp, delTri%list, delTri%lptr, delTri%lend, self%grad3(:,i), errCode)
enddo
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%xComp, delTri%list, delTri%lptr, delTri%lend, self%grad1, sigmaTol, self%sigma1, dSig, errCode)
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%yComp, delTri%list, delTri%lptr, delTri%lend, self%grad2, sigmaTol, self%sigma2, dSig, errCode)
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%zComp, delTri%list, delTri%lptr, delTri%lend, self%grad3, sigmaTol, self%sigma3, dSig, errCode)
end subroutine
subroutine SetScalarSourceData(self, aMesh, delTri, scalarField)
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: scalarField
integer(kint) :: i, errCode
real(kreal) :: dSig
do i = 1, aMesh%particles%N
call GRADL( aMesh%particles%n, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
scalarField%scalar, delTri%list, delTri%lptr, delTri%lend, self%grad1(:,i), errCode)
enddo
call GETSIG(aMesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
scalarField%scalar, delTri%list, delTri%lptr, delTri%lend, self%grad1, sigmaTol, self%sigma1, dSig, errCode)
end subroutine
subroutine SetSourceLagrangianParameter( self, aMesh, delTri )
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
integer(kint) :: i, errCode
real(kreal) :: dSig
do i = 1, aMesh%particles%N
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%x0, delTri%list, delTri%lptr, delTri%lend, self%grad1(:,i), errCode)
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%y0, delTri%list, delTri%lptr, delTri%lend, self%grad2(:,i), errCode)
call GRADL( aMesh%particles%N, i, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%z0, delTri%list, delTri%lptr, delTri%lend, self%grad3(:,i), errCode)
enddo
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%x0, delTri%list, delTri%lptr, delTri%lend, self%grad1, sigmaTol, self%sigma1, dSig, errCode)
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%y0, delTri%list, delTri%lptr, delTri%lend, self%grad2, sigmaTol, self%sigma2, dSig, errCode)
call GETSIG(amesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%z0, delTri%list, delTri%lptr, delTri%lend, self%grad3, sigmaTol, self%sigma3, dSig, errCode)
end subroutine
function InterpolateScalar( lon, lat, self, aMesh, delTri, scalarField)
real(kreal) :: InterpolateScalar
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: scalarField
real(kreal), intent(in) :: lon
real(kreal), intent(in) :: lat
integer(kint) :: errCode
call INTRC1( aMesh%particles%N, lat, lon, amesh%particles%x, amesh%particles%y, amesh%particles%z, &
scalarField%scalar, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma1, GRAD_FLAG, &
self%grad1, startTriangle, InterpolateScalar, errCode)
end function
subroutine InterpolateScalarToUnifLatLonGrid( self, aMesh, delTri, scalarField, lons, lats, interpOut)
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: scalarField
real(kreal), dimension(:), intent(in) :: lons
real(kreal), dimension(:), intent(in) :: lats
real(kreal), dimension(:,:), intent(out) :: interpOut
!
type(MPISetup) :: mpiLons
integer(kint) :: i, j, nLat, nLon, errCode
nLon = size(lons)
nLat = size(lats)
call New(mpiLons, nLon, numProcs)
do j = mpiLons%indexStart(procRank), mpiLons%indexEnd(procrank)
do i = 1, nLat
call INTRC1( aMesh%particles%N, lats(i), lons(j), aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
scalarField%scalar, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma1, GRAD_FLAG, &
self%grad1, startTriangle, interpOut(i,j), errCode)
enddo
enddo
do i = 0, numProcs - 1
call MPI_BCAST( interpOut(:,mpiLons%indexStart(i):mpiLons%indexEnd(i)), nLat * mpiLons%messageLength(i), &
MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, errCode)
enddo
call Delete(mpiLons)
end subroutine
function InterpolateVector( lon, lat, self, aMesh, delTri, vectorField)
real(kreal), dimension(3) :: InterpolateVector
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: vectorField
real(kreal), intent(in) :: lon
real(kreal), intent(in) :: lat
integer(kint) :: errCode
call INTRC1( aMesh%particles%n, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%xComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma1, GRAD_FLAG, &
self%grad1, startTriangle, InterpolateVector(1), errCode)
call INTRC1( aMesh%particles%n, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%yComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma2, GRAD_FLAG, &
self%grad2, startTriangle, InterpolateVector(2), errCode)
call INTRC1( aMesh%particles%n, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%zComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma3, GRAD_FLAG, &
self%grad3, startTriangle, InterpolateVector(3), errCode)
end function
subroutine InterpolateVectorToUnifLatLonGrid( self, aMesh, delTri, vectorField, lons, lats, interpX, interpY, interpZ)
type(SSRFPACKInterface), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
type(Field), intent(in) :: vectorField
real(kreal), dimension(:), intent(in) :: lons
real(kreal), dimension(:), intent(in) :: lats
real(kreal), dimension(:,:), intent(out) :: interpX
real(kreal), dimension(:,:), intent(out) :: interpY
real(kreal), dimension(:,:), intent(out) :: interpZ
!
integer(kint) :: i, j, errCode, nLat, nLon
type(MPISetup) :: mpiLons
nLat = size(lats)
nLon = size(lons)
call New(mpiLons, nLon, numProcs)
do j = mpiLons%indexStart(procRank), mpiLons%indexEnd(procRank)
do i = 1, nLat
call INTRC1( aMesh%particles%n, lats(i), lons(j), aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%xComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma1, GRAD_FLAG, &
self%grad1, startTriangle, interpX(i,j), errCode)
call INTRC1( aMesh%particles%n, lats(i), lons(j), aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%yComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma2, GRAD_FLAG, &
self%grad2, startTriangle, interpY(i,j), errCode)
call INTRC1( aMesh%particles%n, lats(i), lons(j), aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
vectorField%zComp, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma3, GRAD_FLAG, &
self%grad3, startTriangle, interpZ(i,j), errCode)
enddo
enddo
do i = 0, numProcs - 1
call MPI_BCAST( interpX(:,mpiLons%indexstart(i):mpiLons%indexEnd(i)), nLat * mpiLons%messageLength(i), &
MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, errCode)
call MPI_BCAST( interpY(:,mpiLons%indexstart(i):mpiLons%indexEnd(i)), nLat * mpiLons%messageLength(i), &
MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, errCode)
call MPI_BCAST( interpZ(:,mpiLons%indexstart(i):mpiLons%indexEnd(i)), nLat * mpiLons%messageLength(i), &
MPI_DOUBLE_PRECISION, i, MPI_COMM_WORLD, errCode)
enddo
call Delete(mpiLons)
end subroutine
function InterpolateLagParam(lon, lat, self, aMesh, delTri )
real(kreal), dimension(3) :: InterpolateLagParam
type(SSRFPACKInterface), intent(in) :: self
type(PolyMesh2d), intent(in) :: aMesh
type(DelaunayTriangulation), intent(in) :: delTri
real(kreal), intent(in) :: lon
real(kreal), intent(in) :: lat
integer(kint) :: errCode
call INTRC1( aMesh%particles%N, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%x0, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma1, GRAD_FLAG, &
self%grad1, startTriangle, InterpolateLagParam(1), errCode)
call INTRC1( aMesh%particles%N, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%y0, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma2, GRAD_FLAG, &
self%grad2, startTriangle, InterpolateLagParam(2), errCode)
call INTRC1( aMesh%particles%N, lat, lon, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
aMesh%particles%z0, delTri%list, delTri%lptr, delTri%lend, SIGMA_FLAG, self%sigma3, GRAD_FLAG, &
self%grad3, startTriangle, InterpolateLagParam(3), errCode)
end function
!
!----------------
! private methods
!----------------
!
subroutine BuildDelaunayTriangulation(self, aMesh)
type(DelaunayTriangulation), intent(inout) :: self
type(PolyMesh2d), intent(in) :: aMesh
!
real(kreal), allocatable, dimension(:) :: dist
integer(kint), allocatable, dimension(:) :: near
integer(kint), allocatable, dimension(:) :: next
integer(kint) :: lnew, errCode
allocate(dist(amesh%particles%n))
allocate(near(amesh%particles%n))
allocate(next(amesh%particles%n))
call TRMESH(aMesh%particles%N, aMesh%particles%x, aMesh%particles%y, aMesh%particles%z, &
self%list, self%lptr, self%lend, lnew, near, next, dist, errCode )
if ( errCode == -1 ) then
call LogMessage(log,ERROR_LOGGING_LEVEL, trim(logKey)//' TRMESH ERROR :',' found n < 3 points.')
elseif (errCode == -2 ) then
call LogMessage(log,ERROR_LOGGING_LEVEL,trim(logKey)//' TRMESH ERROR :',' found first three nodes to be colinear.')
elseif ( errCode > 0 ) then
write(logString,'(A,I8,A)') ' node ',errCode,' is a duplicate.'
call LogMessage(log,ERROR_LOGGING_LEVEL,trim(logKey)//' TRMESH ERROR :', trim(logString))
endif
deallocate(dist)
deallocate(near)
deallocate(next)
end subroutine
subroutine InitLogger(aLog,rank)
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