-
Notifications
You must be signed in to change notification settings - Fork 1
/
BIVARParticlesInterface.f90
169 lines (137 loc) · 5.89 KB
/
BIVARParticlesInterface.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
module BIVARInterfaceModule
use NumberKindsModule
use LoggerModule
use ParticlesModule
use FieldModule
use BIVARModule
implicit none
private
public BIVARInterface, New, Delete
public InterpolateScalar, InterpolateVector
!
!----------------
! Types and module constants
!----------------
!
type BIVARInterface
real(kreal), dimension(:), pointer :: estPartials1 => null()
real(kreal), dimension(:), pointer :: estPartials2 => null()
end type
type PlaneTri
integer(kint) :: nTri
integer(kint), pointer, dimension(:) :: vertIndices => null()
integer(kint) :: nBorderSegments
integer(kint), pointer, dimension(:) :: borderIndices => null()
integer(kint), pointer, dimension(:) :: integerWork => null()
real(kreal), pointer, dimension(:) :: realWork => null()
contains
final :: deleteTri
end type
interface New
module procedure newPrivateField
module procedure newPrivateLagParam
module procedure newTri
end interface
interface Delete
module procedure deletePrivate
module procedure deleteTri
end interface
contains
subroutine newTri( delTri, aParticles )
type(PlaneTri), intent(out) :: delTri
type(Particles), intent(in) :: aParticles
!
integer(kint), allocatable, dimension(:) :: iwp
allocate(delTri%vertIndices( 6 * aParticles%N - 15) )
allocate(delTri%borderIndices(6* aParticles%N) )
allocate(delTri%integerWork( 18 * aParticles%N))
allocate(delTri%realWork(aParticles%N))
allocate(iwp( aParticles%N))
call IDTANG( aParticles%N, aParticles%x, aParticles%y, delTri%nTri, delTri%vertIndices, delTri%nBorderSegments, &
delTri%borderIndices, delTri%integerWork, iwp, delTri%realWork)
deallocate(iwp)
end subroutine
subroutine deleteTri(delTri)
type(PlaneTri), intent(inout) :: delTri
if ( associated(delTri%vertIndices)) deallocate(delTri%vertIndices)
if ( associated(delTri%borderIndices)) deallocate(delTri%borderIndices)
if ( associated(delTri%integerWork)) deallocate(delTri%integerWork)
if ( associated(delTri%realWork)) deallocate(delTri%realWork)
end subroutine
subroutine newPrivateField( self, delTri, sourceParticles, sourceField )
type(BIVARInterface), intent(out) :: self
type(PlaneTri), intent(in) :: delTri
type(Particles), intent(in) :: sourceParticles
type(Field), intent(in) :: sourceField
!
real(kreal), allocatable, dimension(:) :: realWork
allocate(self%estPartials1( 5 * sourceParticles%N ))
if ( sourceField%nDim == 2 ) allocate(self%estPartials2( 5 * sourceParticles%N))
allocate(realWork(sourceParticles%N))
if ( sourceField%nDim == 1 ) then
call IDPDRV( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%scalar, &
delTri%nTri, delTri%vertIndices, self%estPartials1, realWork)
elseif ( sourceField%nDim == 2 ) then
call IDPDRV( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%xComp, &
delTri%nTri, delTri%vertIndices, self%estPartials1, realWork)
call IDPDRV( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%yComp, &
delTri%nTri, delTri%vertIndices, self%estPartials1, realWork)
endif
deallocate(realWork)
end subroutine
subroutine newPrivateLagParam( self, delTri, sourceParticles )
type(BIVARInterface), intent(out) :: self
type(PlaneTri), intent(in) :: delTri
type(Particles), intent(in) :: sourceParticles
!
real(kreal), allocatable, dimension(:) :: realWork
allocate(self%estPartials1( 5 * sourceParticles%N ))
allocate(self%estPartials2( 5 * sourceParticles%N ))
allocate(realWork(sourceParticles%N))
call IDPDRV( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceParticles%x0, &
delTri%nTri, delTri%vertIndices, self%estPartials1, realWork)
call IDPDRV( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceParticles%y0, &
delTri%nTri, delTri%vertIndices, self%estPartials1, realWork)
deallocate(realWork)
end subroutine
subroutine deletePrivate(self)
type(BIVARInterface), intent(inout) :: self
if ( associated(self%estPartials1)) deallocate(self%estPartials1)
if ( associated(self%estPartials2)) deallocate(self%estPartials2)
end subroutine
function InterpolateScalar( self, delTri, sourceParticles, sourceField, xOut, yOut )
real(kreal) :: InterpolateScalar
type(BIVARInterface), intent(in) :: self
type(PlaneTri), intent(in) :: delTri
type(Particles), intent(in) :: sourceParticles
type(Field), intent(in) :: sourceField
real(kreal), intent(in) :: xOut
real(kreal), intent(in) :: yOut
!
integer(kint) :: inTri
call IDLCTN( sourceParticles%N, sourceParticles%x, sourceParticles%y, delTri%nTri, delTri%vertIndices, &
delTri%nBorderSegments, delTri%borderIndices, xOut, yOut, inTri, delTri%integerWork, delTri%realWork)
call IDPTIP( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%scalar, &
delTri%nTri, delTri%vertIndices, delTri%nBorderSegments, delTri%borderIndices, self%estPartials1, inTri, &
xOut, yOut, InterpolateScalar)
end function
function InterpolateVector( self, delTri, sourceParticles, sourceField, xOut, yOut )
real(kreal), dimension(2) :: InterpolateVector
type(BIVARInterface), intent(in) :: self
type(PlaneTri), intent(in) :: delTri
type(Particles), intent(in) :: sourceParticles
type(Field), intent(in) :: sourceField
real(kreal), intent(in) :: xOut
real(kreal), intent(in) :: yOut
!
integer(kint) :: inTri
call IDLCTN( sourceParticles%N, sourceParticles%x, sourceParticles%y, delTri%nTri, delTri%vertIndices, &
delTri%nBorderSegments, delTri%borderIndices, xOut, yOut, inTri, delTri%integerWork, delTri%realWork)
call IDPTIP( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%xComp, &
delTri%nTri, delTri%vertIndices, delTri%nBorderSegments, delTri%borderIndices, self%estPartials1, inTri, &
xOut, yOut, InterpolateVector(1))
call IDPTIP( sourceParticles%N, sourceParticles%x, sourceParticles%y, sourceField%yComp, &
delTri%nTri, delTri%vertIndices, delTri%nBorderSegments, delTri%borderIndices, self%estPartials2, inTri, &
xOut, yOut, InterpolateVector(2))
end function
end module