-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathclasses.txt
345 lines (267 loc) · 11.8 KB
/
classes.txt
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
module classes
use precision
use MpiUtils
use interpolation, only : TSpline1D, TCubicSpline, TLogRegularCubicSpline
use IniObjects
implicit none
Type MatterTransferData
!Computed data
integer :: num_q_trans ! number of steps in k for transfer calculation
real(dl), dimension (:), allocatable :: q_trans
real(dl), dimension (:), allocatable :: sigma_8
real(dl), dimension (:), allocatable :: sigma2_vdelta_8 !growth from sigma_{v delta}
real, dimension(:,:,:), allocatable :: TransferData
!TransferData(entry,k_index,z_index) for entry=Tranfer_kh.. Transfer_tot
contains
procedure :: Free => MatterTransferData_Free
end Type MatterTransferData
Type MatterPowerData
!everything is a function of k/h
integer :: num_k, num_z
real(dl), dimension(:), allocatable :: log_kh, redshifts
!matpower is log(P_k)
real(dl), dimension(:,:), allocatable :: matpower, ddmat
!if NonLinear, nonlin_ratio = sqrt(P_nonlinear/P_linear)
!function of k and redshift NonLinearScaling(k_index,z_index)
real(dl), dimension(:,:), allocatable :: nonlin_ratio
!Sources
real(dl), dimension(:), allocatable :: log_k
real(dl), dimension(:,:), allocatable :: vvpower, ddvvpower
real(dl), dimension(:,:), allocatable :: vdpower, ddvdpower
real(dl), dimension(:,:), allocatable :: nonlin_ratio_vv
real(dl), dimension(:,:), allocatable :: nonlin_ratio_vd
end Type MatterPowerData
!Classes that can be accessed from python and contain allocatable objects and arrays or other classes
!In Python inherited from F2003Class (defined in baseconfig)
!All python-accessible inherited classes must define SelfPointer, and use @fortran_class decorator in python
!Allocatable objects can only contain instances of python-accessible classes so they can be identified from python.
!Class could be abstract and SelfPointer deferred, but Fortran then doesn't allow called inherited "super()"
!methods implemented in abstract classes
Type TPythonInterfacedClass
contains
!SelfPointer msust be overriden for each class to be referenced from python.
!Gets a class pointer from an untyped pointer.
procedure, nopass :: SelfPointer
!PythonClass gets string of python class name; not actually used internally
procedure, nopass :: PythonClass
!Replace is for copying over this type with an instance of the same type
!It must be defined for any class used as a non-allocatable python structure field that can be assigned to
procedure :: Replace
end type TPythonInterfacedClass
Type PythonClassPtr
class(TPythonInterfacedClass), pointer :: Ref => null()
end type
Type PythonClassAllocatable
class(TPythonInterfacedClass), allocatable :: P
end Type PythonClassAllocatable
type, extends(TPythonInterfacedClass) :: TCambComponent
contains
procedure :: ReadParams => TCambComponent_ReadParams
procedure :: Validate => TCambComponent_Validate
end type TCambComponent
type, extends(TPythonInterfacedClass) :: TCAMBParameters
!Actual type defined in model.f90
end type TCAMBParameters
Type, extends(TPythonInterfacedClass) :: TCAMBdata
!Actual type defined in results.f90
end type TCAMBdata
type, extends(TCambComponent) :: TNonLinearModel
real(dl) :: Min_kh_nonlinear = 0.0_dl
contains
procedure :: Init => TNonLinearModel_init
procedure :: GetNonLinRatios => TNonLinearModel_GetNonLinRatios
procedure :: GetNonLinRatios_All => TNonLinearModel_GetNonLinRatios_All
end type TNonLinearModel
type, extends(TCambComponent) :: TInitialPower
contains
procedure :: ScalarPower => TInitialPower_ScalarPower
procedure :: TensorPower => TInitialPower_TensorPower
procedure :: Init => TInitialPower_Init
procedure :: Effective_ns => TInitalPower_Effective_ns
end type TInitialPower
Type, extends(TCambComponent) :: TRecombinationModel
real(dl) :: min_a_evolve_Tm = 1/(1+900.) !scale factor at which to start evolving Delta_TM
contains
procedure :: Init => TRecombinationModel_init
procedure :: x_e => TRecombinationModel_xe !ionization fraction
procedure :: xe_Tm => TRecombinationModel_xe_Tm !ionization fraction and baryon temperature
procedure :: T_m => TRecombinationModel_tm !baryon temperature
procedure :: T_s => TRecombinationModel_ts !Spin temperature
procedure :: Version => TRecombinationModel_version
procedure :: dDeltaxe_dtau => TRecombinationModel_dDeltaxe_dtau
procedure :: get_Saha_z => TRecombinationModel_Get_Saha_z
end type
Type, extends(TCambComponent) :: TReionizationModel
logical :: Reionization = .true.
contains
procedure :: Init => TReionizationModel_Init
procedure :: x_e => TReionizationModel_xe
procedure :: get_timesteps => TReionizationModel_get_timesteps
end Type TReionizationModel
interface
subroutine TClassDverk (this,n, fcn, x, y, xend, tol, ind, c, nw, w)
use Precision
import
class(TCambComponent), target :: this
integer n, ind
real(dl) x, y(n), xend, tol, c(*), w(nw,9)
external fcn
end subroutine
end interface
contains
subroutine TCambComponent_ReadParams(this, Ini)
class(TCambComponent) :: this
class(TIniFile), intent(in) :: Ini
end subroutine TCambComponent_ReadParams
subroutine TCambComponent_Validate(this, OK)
class(TCambComponent),intent(in) :: this
logical, intent(inout) :: OK
end subroutine TCambComponent_Validate
function PythonClass()
character(LEN=:), allocatable :: PythonClass
PythonClass = ''
error stop 'PythonClass Not implemented'
end function PythonClass
subroutine SelfPointer(cptr, P)
use iso_c_binding
Type(c_ptr) :: cptr
Type (TPythonInterfacedClass), pointer :: PType
class (TPythonInterfacedClass), pointer :: P
call c_f_pointer(cptr, PType)
P => PType
error stop 'Class must define SelfPointer function returning pointer to actual type'
end subroutine SelfPointer
subroutine Replace(this, replace_with)
class(TPythonInterfacedClass), target :: this
class(TPythonInterfacedClass) :: replace_with
error stop 'Assignment not implemented for this class'
end subroutine Replace
subroutine TNonLinearModel_Init(this, State)
class(TNonLinearModel) :: this
class(TCAMBdata), target :: State
end subroutine TNonLinearModel_Init
subroutine TNonLinearModel_GetNonLinRatios(this,State,CAMB_Pk)
class(TNonLinearModel) :: this
class(TCAMBdata) :: State
type(MatterPowerData), target :: CAMB_Pk
error stop 'GetNonLinRatios Not implemented'
end subroutine TNonLinearModel_GetNonLinRatios
subroutine TNonLinearModel_GetNonLinRatios_All(this,State,CAMB_Pk)
class(TNonLinearModel) :: this
class(TCAMBdata) :: State
type(MatterPowerData), target :: CAMB_Pk
error stop 'GetNonLinRatios_all not supported (no non-linear velocities)'
end subroutine TNonLinearModel_GetNonLinRatios_All
function TInitialPower_ScalarPower(this, k)
class(TInitialPower) :: this
real(dl), intent(in) ::k
real(dl) TInitialPower_ScalarPower
TInitialPower_ScalarPower = 0
error stop 'ScalarPower not implemented'
end function TInitialPower_ScalarPower
function TInitialPower_TensorPower(this, k)
class(TInitialPower) :: this
real(dl), intent(in) ::k
real(dl) TInitialPower_TensorPower
TInitialPower_TensorPower = 0
error stop 'TensorPower not implemented'
end function TInitialPower_TensorPower
subroutine TInitialPower_Init(this, Params)
class(TInitialPower) :: this
class(TCAMBParameters), intent(in) :: Params
end subroutine TInitialPower_Init
function TInitalPower_Effective_ns(this)
class(TInitialPower) :: this
real(dl) :: TInitalPower_Effective_ns
TInitalPower_Effective_ns=1
call MpiStop('Power spectrum does not support an effective n_s for halofit')
end function TInitalPower_Effective_ns
subroutine MatterTransferData_Free(this)
class(MatterTransferData):: this
integer st
deallocate(this%q_trans, STAT = st)
deallocate(this%TransferData, STAT = st)
deallocate(this%sigma_8, STAT = st)
deallocate(this%sigma2_vdelta_8, STAT = st)
end subroutine MatterTransferData_Free
function TRecombinationModel_tm(this,a)
class(TRecombinationModel) :: this
real(dl), intent(in) :: a
real(dl) TRecombinationModel_tm
call MpiStop('TRecombinationModel_tm not implemented')
TRecombinationModel_tm=0
end function TRecombinationModel_tm
function TRecombinationModel_ts(this,a)
class(TRecombinationModel) :: this
!zrec(1) is zinitial-delta_z
real(dl), intent(in) :: a
real(dl) TRecombinationModel_ts
call MpiStop('TRecombinationModel_ts not implemented')
TRecombinationModel_ts=0
end function TRecombinationModel_ts
function TRecombinationModel_xe(this,a)
class(TRecombinationModel) :: this
real(dl), intent(in) :: a
real(dl) TRecombinationModel_xe
call MpiStop('TRecombinationModel_xe not implemented')
TRecombinationModel_xe=0
end function TRecombinationModel_xe
subroutine TRecombinationModel_xe_Tm(this,a, xe, Tm)
!Not required to implement, but may be able to optimize
class(TRecombinationModel) :: this
real(dl), intent(in) :: a
real(dl), intent(out) :: xe, Tm
xe = this%x_e(a)
Tm = this%T_m(a)
end subroutine TRecombinationModel_xe_Tm
function TRecombinationModel_version(this) result(version)
class(TRecombinationModel) :: this
character(LEN=:), allocatable :: version
version = ''
end function TRecombinationModel_version
subroutine TRecombinationModel_init(this,State, WantTSpin)
class(TRecombinationModel), target :: this
class(TCAMBdata), target :: State
logical, intent(in), optional :: WantTSpin
end subroutine TRecombinationModel_init
function TRecombinationModel_dDeltaxe_dtau(this,a, Delta_xe,Delta_nH, Delta_Tm, hdot, kvb,adotoa)
!d x_e/d tau
class(TRecombinationModel) :: this
real(dl) TRecombinationModel_dDeltaxe_dtau
real(dl), intent(in):: a, Delta_xe,Delta_nH, Delta_Tm, hdot, kvb,adotoa
call MpiStop('TRecombinationModel_dDeltaxe_dtau not implemented')
TRecombinationModel_dDeltaxe_dtau=0
end function TRecombinationModel_dDeltaxe_dtau
real(dl) function TRecombinationModel_Get_Saha_z(this)
class(TRecombinationModel) :: this
call MpiStop('TRecombinationModel_Get_Saha_z not implemented')
TRecombinationModel_Get_Saha_z = 0
end function
function TReionizationModel_xe(this, z, tau, xe_recomb)
!a and time tau and redundant, both provided for convenience
!xe_recomb is xe(tau_start) from recombination (typically very small, ~2e-4)
!xe should map smoothly onto xe_recomb
class(TReionizationModel) :: this
real(dl), intent(in) :: z
real(dl), intent(in), optional :: tau, xe_recomb
real(dl) TReionizationModel_xe
call MpiStop('TReionizationModel_xe not implemented')
TReionizationModel_xe=0
end function TReionizationModel_xe
subroutine TReionizationModel_get_timesteps(this, n_steps, z_start, z_complete)
!minimum number of time steps to use between tau_start and tau_complete
!Scaled by AccuracyBoost later
!steps may be set smaller than this anyway
class(TReionizationModel) :: this
integer, intent(out) :: n_steps
real(dl), intent(out):: z_start, z_Complete
call MpiStop('TReionizationModel_get_timesteps not implemented')
n_steps=0
z_start=0
z_complete=0
end subroutine TReionizationModel_get_timesteps
subroutine TReionizationModel_Init(this, State)
class(TReionizatioNModel) :: this
class(TCAMBdata), target :: State
end subroutine TReionizationModel_Init
end module classes