-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathbhe_models.py
1917 lines (1642 loc) · 80.8 KB
/
bhe_models.py
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
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
from __future__ import division
import numpy as np
import math
from scipy import sparse
from scipy.sparse import csr_matrix
from scipy.optimize import minimize_scalar
class BHE_CXA_impl:
'''
implicit simulation model for coaxial bhe with anular inlet
'''
def setTimestep(self,dt):
self.dt = dt
def setnz(self,nz):
self.nz = nz
def setSoilBC(self,Tsoil):
self.result[3*self.nz:4*self.nz] = Tsoil
def getGroutBC(self):
return np.mean(self.T_grout[1:])
def getFluidOut(self):
return self.Tf_out[1]
def getFluidIn(self):
return self.Tf_in[1]
def initialize(self,BheData):
self.dynviscF = BheData['dynviscF']
self.dz = BheData['length']/self.nz
self.nz = self.nz+1
# Geometry
self.ro = BheData['diamB']/2
self.ri_i = BheData['odiamPin']/2-BheData['thickPin']
self.ri_o = BheData['odiamPin']/2
self.ro_i = BheData['odiamPout']/2-BheData['thickPout']
self.ro_o = BheData['odiamPout']/2
self.di_i = BheData['odiamPin']-2*BheData['thickPin']
self.do_i = BheData['odiamPout']-2*BheData['thickPout']
self.dh = self.di_i-BheData['odiamPout']
# Flow velocity
self.uo = BheData['Qf']/(np.pi*self.ro_i**2)
self.ui = BheData['Qf']/(np.pi*(self.ri_i**2-self.ro_o**2))
self.maxVel = np.max([self.ui,self.uo])
# Flow parameters
self.Pr = self.dynviscF*BheData['capF']/BheData['densF']/BheData['lmF']
self.Reo = self.uo*self.do_i/(self.dynviscF/BheData['densF'])
self.Rei = self.ui*self.dh/(self.dynviscF/BheData['densF'])
self.Nuo = NusseltCoaxo (self.Reo,self.Pr,self.ro_i,BheData['length'])
self.Nui = NusseltCoaxi (self.Rei,self.Pr,BheData['odiamPout'],self.di_i,self.dh,BheData['length'])
self.Radvo = 1/(self.Nuo*BheData['lmF']*np.pi)
self.Radvia = 1/(self.Nui*BheData['lmF']*np.pi) * self.dh/BheData['odiamPout']
self.Radvib = 1/(self.Nui*BheData['lmF']*np.pi) * self.dh/self.di_i
# Resistances
self.x = np.log(((BheData['diamB']**2+BheData['odiamPin']**2)**0.5)/((2**0.5)*BheData['odiamPin']))/np.log(BheData['diamB']/BheData['odiamPin'])
self.Rg = np.log(BheData['diamB']/BheData['odiamPin'])/(2*np.pi*BheData['lmG'])
self.Rgs = (1-self.x)*self.Rg
self.Rgs_coupling = 1./self.Rgs
self.Rconb = self.x*self.Rg
self.Rcono = np.log(self.ro_o/self.ro_i)/(2*np.pi*BheData['lmPout'])
self.Rconi = np.log(self.ri_o/self.ri_i)/(2*np.pi*BheData['lmPin'])
self.Rff = self.Radvo + self.Radvia + self.Rcono
self.Rfig = self.Radvib + self.Rconi + self.Rconb
self.RffNoFlow = self.Rcono
self.RfigNoFlow = self.Rconb + self.Rconi
# Volumes and Areas
self.Ag = np.pi*(self.ro**2-self.ri_o**2) # area grout
self.Vg = self.Ag*self.dz # volume grout
self.Ain = np.pi*(self.ri_i**2-self.ro_o**2) # area pipeIn/fluidIn
self.Vin = self.Ain*self.dz # volume pipeIn/fluidIn
self.Aout = np.pi*self.ro_i**2 # area pipeOut/fluidOut
self.Vout = self.Aout*self.dz # volume pipeOut/fluidOut
# Variables Flow
self.F1 = self.dz/self.Rgs
self.F2 = self.dz/self.Rfig
self.F3 = self.dt/self.Vg/BheData['capG']
self.F4 = self.dz/self.Rff
self.F5i = self.dt/self.dz*self.ui
self.F5o = self.dt/self.dz*self.uo
self.F6 = self.dt/self.Vin/BheData['capF']
self.F7 = self.dt/self.Vout/BheData['capF']
self.F8 = BheData['lmG']/self.dz*self.Ag
self.F6F2 = self.F6*self.F2
self.F6F3 = self.F6*self.F3
self.F6F4 = self.F6*self.F4
self.F7F4 = self.F7*self.F4
self.F3F1 = self.F3*self.F1
self.F3F2 = self.F3*self.F2
self.F3F8 = self.F3*self.F8
# Variables NoFlow
self.F2NoFlow = self.dz/self.RfigNoFlow
self.F4NoFlow = self.dz/self.RffNoFlow
self.F6F2NoFlow = self.F6*self.F2NoFlow
self.F6F4NoFlow = self.F6*self.F4NoFlow
self.F7F4NoFlow = self.F7*self.F4NoFlow
self.F3F2NoFlow = self.F3*self.F2NoFlow
self.F5iNoFlow = 0
self.F5oNoFlow = 0
# Set up matrix Flow
# data = data array, posX = X Position, posY = Y Position
self.na = 4 # number of cell arrays, Tin, Tout, Tgrout, Tsoil
##### Grout
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2)
self.dataTemp[0] = 1+self.F3F8+self.F3F1+self.F3F2 # oberste Zelle
self.dataTemp[self.nz-1] = 1+self.F3F8+self.F3F1+self.F3F2 # unterste Zelle
self.data = self.dataTemp
self.posX = np.arange(0,self.nz,1)
self.posY = np.arange(0,self.nz,1)
## grout z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1))
## grout z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1))
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1)
## grout to soil
self.dataTemp = np.ones(self.nz)*-self.F3F1
self.dataTemp[0] = 0
self.data = np.append(self.data, self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-1)*self.nz)
## grout to tin
self.dataTemp = np.ones(self.nz)*-self.F3F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## main diagonal Tin
self.dataTemp = np.ones(self.nz)*(1+self.F5i+self.F6F2+self.F6F4)
#self.dataTemp[0] = (1+self.F6F2+self.F6F4)
self.dataTemp[0] = 1
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + self.nz)
## Tin z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5i)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+ self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+ self.nz)
## Tin to grout
self.dataTemp = np.ones(self.nz)*-self.F6F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
## Tin to Tout
self.dataTemp = np.ones(self.nz)*-self.F6F4
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-2)*self.nz)
## main diagonal Tout
self.data = np.append(self.data,np.ones(self.nz)*(1+self.F5o+self.F7F4))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 2*self.nz)
## Tout z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5o)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1+2*self.nz)
## Tout to tin
self.dataTemp = np.ones(self.nz)*-self.F7F4
self.dataTemp[self.nz-1] = -(self.F5o+self.F7F4)
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## main diagonal Soil
self.data = np.append(self.data,np.ones(self.nz))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 3*self.nz)
self.K_Flow = csr_matrix((self.data, (self.posX, self.posY)), shape=(self.nz*self.na, self.nz*self.na),dtype=np.float)
self.K_sparse_Flow = sparse.csr_matrix(self.K_Flow)
# Set up matrix NoFlow
# data = Datenarray, posX = X Position, posY = Y Position, dataTemp = Temporär
self.na = 4 # number of cell arrays, Tin, Tout, Tgrout, Tsoil
##### Grout
## main diagonal grout
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2NoFlow)
self.dataTemp[0] = 1+self.F3F8+self.F3F1+self.F3F2NoFlow # oberste Zelle
self.dataTemp[self.nz-1] = 1+self.F3F8+self.F3F1+self.F3F2NoFlow # unterste Zelle
self.data = self.dataTemp
self.posX = np.arange(0,self.nz,1)
self.posY = np.arange(0,self.nz,1)
## grout z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1))
## grout z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1))
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1)
## grout to soil
self.data = np.append(self.data,np.ones(self.nz)*-self.F3F1)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-1)*self.nz)
## grout to tin
self.data = np.append(self.data,np.ones(self.nz)*-self.F3F2NoFlow)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## main diagonal Tin
self.dataTemp = np.ones(self.nz)*(1+self.F5iNoFlow+self.F6F2NoFlow+self.F6F4NoFlow)
#self.dataTemp[0] = (1+self.F6F2NoFlow+self.F6F4NoFlow)
self.dataTemp[0] = 1
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + self.nz)
## Tin z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5iNoFlow)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+ self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+ self.nz)
## Tin to grout
self.data = np.append(self.data,np.ones(self.nz)*-self.F6F2NoFlow)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
## Tin to Tout
self.data = np.append(self.data,np.ones(self.nz)*-self.F6F4NoFlow)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-2)*self.nz)
## main diagonal Tout
self.data = np.append(self.data,np.ones(self.nz)*(1+self.F5oNoFlow+self.F7F4NoFlow))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 2*self.nz)
## Tout z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5oNoFlow)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1+2*self.nz)
## Tout to tin
self.dataTemp = np.ones(self.nz)*-self.F7F4NoFlow
self.dataTemp[self.nz-1] = -(self.F5oNoFlow+self.F7F4NoFlow)
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## main diagonal Soil
self.data = np.append(self.data,np.ones(self.nz))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 3*self.nz)
self.K_NoFlow = csr_matrix((self.data, (self.posX, self.posY)), shape=(self.nz*self.na, self.nz*self.na),dtype=np.float)
self.K_sparse_NoFlow = sparse.csr_matrix(self.K_NoFlow)
######## Told ########
######################
self.result = np.ones(self.nz*self.na)*BheData['Tundist'] # Initialbedingung
self.Tf_in = np.full(self.nz,BheData['Tundist']) # CellArray Tf_in
self.Tf_out = np.full(self.nz,BheData['Tundist']) # CellArray Tf_out
self.T_grout = np.full(self.nz,BheData['Tundist']) # CellArray Tf_grout
return True
def calcSondeFlow(self,tfinal,T_in):
for i in range(0,tfinal):
self.result[self.nz] = T_in # BC Tin
self.U =sparse.linalg.spsolve(self.K_sparse_Flow,self.result) # Solve EQS
# Rausschreiben für Plot
self.T_grout = self.U[0:self.nz]
self.Tf_in = self.U[self.nz:2*self.nz]
self.Tf_out = self.U[2*self.nz:3*self.nz]
# reinschreiben nächster Zeitschritt
self.result[:] = np.copy(self.U[:])
return True
def calcSondeNoFlow(self,tfinal,T_in):
for i in range(0,tfinal):
# Berechnung grout
self.result[self.nz] = T_in # Vorgabe Randbedingung Tin
self.U =sparse.linalg.spsolve(self.K_sparse_NoFlow,self.result) # Lösen GlS
# Rausschreiben für Plot
self.T_grout = self.U[0:self.nz]
self.Tf_in = self.U[self.nz:2*self.nz]
self.Tf_out = self.U[2*self.nz:3*self.nz]
# reinschreiben nächster Zeitschritt
self.result[:] = np.copy(self.U[:])
return True
class BHE_CXA_expl:
'''
explicit simulation model for coaxial bhe with anular inlet
'''
def setTimestep(self,dt):
self.dt = dt
def setnz(self,nz):
self.nz = nz
def setSoilBC(self,Tsoil):
self.TsA[:] = Tsoil
def getGroutBC(self):
return np.mean(self.T_grout)
def getFluidOut(self):
return self.Tf_out[0]
def getFluidIn(self):
return self.Tf_in[0]
def initialize(self,BheData):
self.dz = BheData['length']/self.nz
dynviscF = BheData['dynviscF']
# Geometry
D = BheData['diamB']
self.ro = BheData['diamB']/2
self.ri_i = BheData['odiamPin']/2-BheData['thickPin']
self.ri_o = BheData['odiamPin']/2
self.ro_i = BheData['odiamPout']/2-BheData['thickPout']
self.ro_o = BheData['odiamPout']/2
self.di_i = BheData['odiamPin']-2*BheData['thickPin']
self.do_i = BheData['odiamPout']-2*BheData['thickPout']
self.dh = self.di_i-BheData['odiamPout']
# Flow velocity
self.uo = BheData['Qf']/(np.pi*self.ro_i**2)
self.ui = BheData['Qf']/(np.pi*(self.ri_i**2-self.ro_o**2))
self.maxVel = np.max([self.ui,self.uo])
# Flow parameters
self.Pr = BheData['dynviscF']*BheData['capF']/BheData['densF']/BheData['lmF']
self.Reo = self.uo*self.do_i/(BheData['dynviscF']/BheData['densF'])
self.Rei = self.ui*self.dh/(BheData['dynviscF']/BheData['densF'])
self.Nuo = NusseltCoaxo (self.Reo,self.Pr,self.ro_i,BheData['length'])
self.Nui = NusseltCoaxi (self.Rei,self.Pr,BheData['odiamPout'],self.di_i,self.dh,BheData['length'])
self.Radvo = 1/(self.Nuo*BheData['lmF']*np.pi)
self.Radvia = 1/(self.Nui*BheData['lmF']*np.pi) * self.dh/BheData['odiamPout']
self.Radvib = 1/(self.Nui*BheData['lmF']*np.pi) * self.dh/self.di_i
# Resistances
self.x = np.log(((BheData['diamB']**2+BheData['odiamPin']**2)**0.5)/((2**0.5)*BheData['odiamPin']))/np.log(BheData['diamB']/BheData['odiamPin'])
self.Rg = np.log(BheData['diamB']/BheData['odiamPin'])/(2*np.pi*BheData['lmG'])
self.Rgs = (1-self.x)*self.Rg
self.Rgs_coupling = 1./self.Rgs
self.Rconb = self.x*self.Rg
self.Rcono = np.log(self.ro_o/self.ro_i)/(2*np.pi*BheData['lmPout'])
self.Rconi = np.log(self.ri_o/self.ri_i)/(2*np.pi*BheData['lmPin'])
self.Rff = self.Radvo + self.Radvia + self.Rcono
self.Rfig = self.Radvib + self.Rconi + self.Rconb
self.RffNoFlow = self.Rcono
self.RfigNoFlow = self.Rconb + self.Rconi
# Volumes and Areas
self.Ag = np.pi*(self.ro**2-self.ri_o**2) # area grout
self.Vg = self.Ag*self.dz # volume grout
self.Ain = np.pi*(self.ri_i**2-self.ro_o**2) # area pipeIn/fluidIn
self.Vin = self.Ain*self.dz # volume pipeIn/fluidIn
self.Aout = np.pi*self.ro_i**2 # area pipeOut/fluidOut
self.Vout = self.Aout*self.dz # volume pipeOut/fluidOut
# Variables Flow
self.F1 = self.dz/self.Rgs
self.F2 = self.dz/self.Rfig
self.F3 = self.dt/self.Vg/BheData['capG']
self.F4 = self.dz/self.Rff
self.F5i = self.dt/self.dz*self.ui
self.F5o = self.dt/self.dz*self.uo
self.F6 = self.dt/self.Vin/BheData['capF']
self.F7 = self.dt/self.Vout/BheData['capF']
self.F8 = BheData['lmG']/self.dz*self.Ag
self.F6F2 = self.F6*self.F2
self.F6F3 = self.F6*self.F3
self.F6F4 = self.F6*self.F4
self.F7F4 = self.F7*self.F4
self.F3F1 = self.F3*self.F1
self.F3F2 = self.F3*self.F2
self.F3F8 = self.F3*self.F8
# Variables NoFlow
self.F2NoFlow = self.dz/self.RfigNoFlow
self.F4NoFlow = self.dz/self.RffNoFlow
self.F6F2NoFlow = self.F6*self.F2NoFlow
self.F6F4NoFlow = self.F6*self.F4NoFlow
self.F7F4NoFlow = self.F7*self.F4NoFlow
self.F3F2NoFlow = self.F3*self.F2NoFlow
self.idxPlus = np.roll(np.linspace(0,self.nz-1,self.nz,dtype = 'int'),1)
self.idxMinus = np.roll(np.linspace(0,self.nz-1,self.nz,dtype = 'int'),-1)
self.Q_cond = np.zeros(self.nz)
self.Phi_w = np.zeros(self.nz)
self.Phi_e = np.zeros(self.nz)
# Calc Arrays
self.Tf_in = np.full(self.nz,BheData['Tundist']) # CellArray Tf_in
self.Tf_out = np.full(self.nz,BheData['Tundist']) # CellArray Tf_out
self.T_grout = np.full(self.nz,BheData['Tundist']) # CellArray Tf_grout
self.TsA = np.full(self.nz,BheData['Tundist']) # CellArray Average Soil Temperature
self.zgrid = np.linspace(0,BheData['length'],self.nz)
self.alt_Tf_in = np.full(self.nz,BheData['Tundist']) # CellArray Tf_in
self.alt_Tf_out = np.full(self.nz,BheData['Tundist']) # CellArray Tf_out
self.alt_T_grout = np.full(self.nz,BheData['Tundist']) # CellArray Tf_grout
return True
def calcSondeFlow(self,tfinal,T_in):
for i in range(0,tfinal):
# Berechnung grout
self.Q_cond = (np.take(self.alt_T_grout,self.idxMinus)-2*self.alt_T_grout+np.take(self.alt_T_grout,self.idxPlus))
self.Q_cond[0] = (self.alt_T_grout[1]-self.alt_T_grout[0])
self.Q_cond[self.nz-1] = (self.alt_T_grout[self.nz-2]-self.alt_T_grout[self.nz-1])
self.T_grout = self.alt_T_grout + self.F3F8 * self.Q_cond + self.F3F1*(self.TsA-self.alt_T_grout) + self.F3F2*(self.alt_Tf_in-self.alt_T_grout)
# Berechnung fluid in
self.Phi_w = np.take(self.alt_Tf_in,self.idxPlus)
self.Phi_w[0] = T_in
self.Tf_in = self.alt_Tf_in + self.F5i*(self.Phi_w - self.alt_Tf_in) + self.F6F2*(self.alt_T_grout-self.alt_Tf_in) + self.F6F4*(self.alt_Tf_out-self.alt_Tf_in)
# Berechnung fluid out
self.Phi_w = np.take(self.alt_Tf_out,self.idxMinus)
self.Phi_w[self.nz-1] = self.alt_Tf_in[self.nz-1]
self.Tf_out = self.alt_Tf_out + self.F5o * (self.Phi_w - self.alt_Tf_out) + self.F7F4*(self.alt_Tf_in-self.alt_Tf_out)
# Rückschreiben
self.alt_T_grout[:] = self.T_grout[:]
self.alt_Tf_in[:] = self.Tf_in[:]
self.alt_Tf_out[:] = self.Tf_out[:]
return True
def calcSondeNoFlow(self,tfinal,T_in):
for i in range(0,tfinal):
# Berechnung grout
self.Q_cond = (np.take(self.alt_T_grout,self.idxMinus)-2*self.alt_T_grout+np.take(self.alt_T_grout,self.idxPlus))
self.Q_cond[0] = (self.alt_T_grout[1]-self.alt_T_grout[0])
self.Q_cond[self.nz-1] = (self.alt_T_grout[self.nz-2]-self.alt_T_grout[self.nz-1])
self.T_grout = self.alt_T_grout + self.F3F8 * self.Q_cond + self.F3F1*(self.TsA-self.alt_T_grout) + self.F3F2NoFlow*(self.alt_Tf_in-self.alt_T_grout)
# Berechnung fluid in
self.Tf_in = self.alt_Tf_in + self.F6F2NoFlow*(self.alt_T_grout-self.alt_Tf_in) + self.F6F4NoFlow*(self.alt_Tf_out-self.alt_Tf_in)
# Berechnung fluid out
self.Tf_out = self.alt_Tf_out + self.F7F4NoFlow*(self.alt_Tf_in-self.alt_Tf_out)
# Rückschreiben
self.alt_T_grout[:] = self.T_grout[:]
self.alt_Tf_in[:] = self.Tf_in[:]
self.alt_Tf_out[:] = self.Tf_out[:]
return True
class BHE_2U_impl:
'''
implicit simulation model for 2U bhe
'''
def setSoilBC(self,Tsoil):
self.result[4*self.nz:5*self.nz] = Tsoil
def getGroutBC(self):
return (np.mean(self.T_grout_in[1:])+np.mean(self.T_grout_out[1:]))/2
def getFluidOut(self):
return self.Tf_out[1]
def getFluidIn(self):
return self.Tf_in[1]
def setTimestep(self,dt):
self.dt = dt
def setnz(self,nz):
self.nz = nz
def initialize(self,BheData):
# discretize
self.dz = BheData['length']/self.nz
self.nz = self.nz+1
# Geometry
self.D = BheData['diamB']
self.ro = BheData['diamB']/2
self.rpi = BheData['odiamP']/2-BheData['thickP']
self.rpo = BheData['odiamP']/2
self.dpi = BheData['odiamP']-2*BheData['thickP']
self.dpo = BheData['odiamP']
# fluid properties
self.dynviscF = BheData['dynviscF']
self.densF = BheData['densF']
self.lmF = BheData['lmF']
self.length = BheData['length']
self.capF = BheData['capF']
# Flow velocity
self.u = BheData['Qf']/(2*np.pi*self.rpi**2)
self.Qold = BheData['Qf']
self.maxVel = self.u
# Flow parameters
self.Pr = self.dynviscF*BheData['capF']/BheData['densF']/BheData['lmF']
self.Re = self.u*self.dpi/(self.dynviscF/BheData['densF'])
self.Nu = Nusselt2U (self.Re,self.Pr,self.rpi,BheData['length'])
self.Radv = 1/(self.Nu*BheData['lmF']*np.pi)
# Thermal Resistances
self.s = BheData['distP']*2**0.5
self.x = np.log((self.D**2+4*self.dpo**2)**0.5/(2*2**0.5*self.dpo))/np.log(self.D/(2*self.dpo))
self.Rg = (3.098-4.432*self.s/self.D + 2.364*self.s**2/(self.D**2)) * np.arccosh((self.D**2+self.dpo**2-self.s**2)/(2*self.D*self.dpo))/(2*np.pi*BheData['lmG'])
self.Rgs = (1-self.x)*self.Rg
self.Rconb = self.x*self.Rg
self.Rcona = np.log(self.rpo/self.rpi)/(2*np.pi*BheData['lmP'])
self.Rar1 = np.arccosh((self.s**2-self.dpo**2)/self.dpo**2)/(2*np.pi*BheData['lmG'])
self.Rar2 = np.arccosh((2*self.s**2-self.dpo**2)/self.dpo**2)/(2*np.pi*BheData['lmG'])
self.Rgg1 = (2*self.Rgs*(self.Rar1-2*self.x*self.Rg))/(2*self.Rgs-self.Rar1+2*self.x*self.Rg)
self.Rgg2 = (2*self.Rgs*(self.Rar2-2*self.x*self.Rg))/(2*self.Rgs-self.Rar2+2*self.x*self.Rg)
### Check negative thermal resistances ###
if (((1/self.Rgg1 + 1/(2*self.Rgs))**(-1) <= 0) or ((1/self.Rgg2 + 1/(2*self.Rgs))**(-1) <= 0)):
self.x = 2/3*self.x
self.Rconb = self.x*self.Rg
self.Rgs = (1-self.x)*self.Rg
self.Rgg1 = (2*self.Rgs*(self.Rar1-2*self.x*self.Rg))/(2*self.Rgs-self.Rar1+2*self.x*self.Rg)
self.Rgg2 = (2*self.Rgs*(self.Rar2-2*self.x*self.Rg))/(2*self.Rgs-self.Rar2+2*self.x*self.Rg)
if (((1/self.Rgg1 + 1/(2*self.Rgs))**(-1) <= 0) or ((1/self.Rgg2 + 1/(2*self.Rgs))**(-1) <= 0)):
self.x = self.x * 1/3
self.Rconb = self.x*self.Rg
self.Rgs = (1-self.x)*self.Rg
self.Rgg1 = (2*self.Rgs*(self.Rar1-2*self.x*self.Rg))/(2*self.Rgs-self.Rar1+2*self.x*self.Rg)
self.Rgg2 = (2*self.Rgs*(self.Rar2-2*self.x*self.Rg))/(2*self.Rgs-self.Rar2+2*self.x*self.Rg)
if (((1/self.Rgg1 + 1/(2*self.Rgs))**(-1) <= 0) or ((1/self.Rgg2 + 1/(2*self.Rgs))**(-1) <= 0)):
self.x = 0
self.Rconb = self.x*self.Rg
self.Rgs = (1-self.x)*self.Rg
self.Rgg1 = (2*self.Rgs*(self.Rar1-2*self.x*self.Rg))/(2*self.Rgs-self.Rar1+2*self.x*self.Rg)
self.Rgg2 = (2*self.Rgs*(self.Rar2-2*self.x*self.Rg))/(2*self.Rgs-self.Rar2+2*self.x*self.Rg)
self.Rgs_coupling = 4./self.Rgs
self.Rfg = self.Radv + self.Rcona + self.Rconb
self.RfgNoFlow = self.Rconb + self.Rcona
# Volumes and Areas
self.Ag = np.pi*(self.ro**2-4*self.rpo**2)/4 # area grout
self.Vg = self.Ag*self.dz # volume grout
self.Ap = np.pi*self.rpi**2 # area pipe/fluid
self.Vp = self.Ap*self.dz # volume pipe/fluid
# Variables Flow # Ab hier werden werte zusammengefasst damit rechnung schneller läuft
self.F1 = self.dz/self.Rgs # am besten zurücksubstituieren für verständnis
self.F2 = self.dz/self.Rfg
self.F3 = self.dt/self.Vg/BheData['capG']
self.F4 = self.dz/self.Rgg1
self.F5 = self.dt/self.dz*self.u
self.F6 = self.dt/self.Vp/BheData['capF']
self.F8 = BheData['lmG']/self.dz*self.Ag
self.F9 = 4./self.Rgs
self.F6F2 = self.F6*self.F2
self.F3F1 = self.F3*self.F1
self.F3F2 = self.F3*self.F2
self.F3F4 = self.F3*self.F4
self.F3F8 = self.F3*self.F8
# Variables NoFlow
self.F5NoFlow = 0
self.F2NoFlow = self.dz/self.RfgNoFlow
self.F3F2NoFlow = self.F3*self.F2NoFlow
self.F6F2NoFlow = self.F6*self.F2NoFlow
# flow dependent variabls
self.flow_vari = np.zeros(6)
self.flow_vari[0] = (1 + self.F5 + self.F6F2)
self.flow_vari[1] = (- self.F5)
self.flow_vari[2] = (- self.F6F2)
self.flow_vari[3] = (1 + 2*self.F3F8 + self.F3F1 + self.F3F2 + 2*self.F3F4)
self.flow_vari[4] = (1 + self.F3F8 + self.F3F1 + self.F3F2 + 2*self.F3F4)
self.flow_vari[5] = (- self.F3F2)
# Set up Matrix
self.na = 5 # number of cell arrays, Tin, Tout, Tgi, Tgo, Tsoil
############### Matrix for flow #################
#################################################
##### Tgi
## main diagonal Tgi
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2+2*self.F3F4)
self.dataTemp[0] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # oberste Zelle
self.dataTemp[self.nz-1] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # unterste Zelle
self.data = self.dataTemp
self.posX = np.arange(0,self.nz,1)
self.posY = np.arange(0,self.nz,1)
## Tgi z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1))
## Tgi z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1))
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1)
## Tgi to soil
self.dataTemp = np.ones(self.nz)*-self.F3F1
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-1)*self.nz)
## Tgi to tin
self.dataTemp = np.ones(self.nz)*-self.F3F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## Tgi to Tgo
self.dataTemp = np.ones(self.nz)*-2*self.F3F4
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-4)*self.nz)
##### Tgo
## main diagonal Tgo
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2+2*self.F3F4)
self.dataTemp[0] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # oberste Zelle
self.dataTemp[self.nz-1] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # unterste Zelle
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+self.nz)
## Tgo z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+self.nz)
## Tgo z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1+self.nz)
## Tgo to soil
self.dataTemp = np.ones(self.nz)*-self.F3F1
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+4*self.nz)
## Tgo to tout
self.dataTemp = np.ones(self.nz)*-self.F3F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+3*self.nz)
## Tgo to Tgi
self.dataTemp = np.ones(self.nz)*-2*self.F3F4
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
##### Tfi
## main diagonal Tfi
self.dataTemp = np.ones(self.nz)*(1+self.F5+self.F6F2)
self.dataTemp[0] = 1
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 2*self.nz)
## Tfi z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+ 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+ 2*self.nz)
## Tfi to Tgi
self.dataTemp = np.ones(self.nz)*-self.F6F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+ 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
##### Tfo
## main diagonal Tfo
self.dataTemp = np.ones(self.nz)*(1+self.F5+self.F6F2)
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 3*self.nz)
## Tfo z+1
self.data = np.append(self.data,np.ones(self.nz)*-self.F5)
self.posXTemp = np.arange(0,self.nz,1)+ 3*self.nz
self.posYTemp = np.arange(0,self.nz,1)+1+3*self.nz
self.posYTemp[self.nz-1] = self.posYTemp[self.nz-1]-(self.nz+1)
self.posX = np.append(self.posX,self.posXTemp)
self.posY = np.append(self.posY,self.posYTemp)
## Tfo to Tgo
self.dataTemp = np.ones(self.nz)*-self.F6F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+ 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+self.nz)
# main diagonal Soil 1
self.data = np.append(self.data,np.ones(self.nz))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 4*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 4*self.nz)
self.K = csr_matrix((self.data, (self.posX, self.posY)), shape=(self.nz*self.na, self.nz*self.na),dtype=np.float)
self.K_sparse_Flow = sparse.csr_matrix(self.K)
self.K_arr = self.K_sparse_Flow.toarray()
# get indices of flow dependent variabels in martix
self.vari_indices = []
for vari in self.flow_vari:
self.indx = [self.K_arr == vari]
self.vari_indices.append(self.indx)
############### Matrix for no flow #################
####################################################
self.F5 = 0
#self.F2 = self.dz/self.RfgNoFlow
self.F3F2 = self.F3*self.F2
self.F6F2 = self.F6*self.F2
##### Tgi
## main diagonal Tgi
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2+2*self.F3F4)
self.dataTemp[0] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # oberste Zelle
self.dataTemp[self.nz-1] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # unterste Zelle
self.data = self.dataTemp
self.posX = np.arange(0,self.nz,1)
self.posY = np.arange(0,self.nz,1)
## Tgi z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1))
## Tgi z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1))
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1)
## Tgi to soil
self.dataTemp = np.ones(self.nz)*-self.F3F1
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-1)*self.nz)
## Tgi to tin
self.dataTemp = np.ones(self.nz)*-self.F3F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-3)*self.nz)
## Tgi to Tgo
self.dataTemp = np.ones(self.nz)*-2*self.F3F4
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1))
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+(self.na-4)*self.nz)
##### Tgo
## main diagonal Tgo
self.dataTemp = np.ones(self.nz)*(1+2*self.F3F8+self.F3F1+self.F3F2+2*self.F3F4)
self.dataTemp[0] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # oberste Zelle
self.dataTemp[self.nz-1] = (1+self.F3F8+self.F3F1+self.F3F2+2*self.F3F4) # unterste Zelle
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+self.nz)
## Tgo z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+self.nz)
## Tgo z+1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F3F8)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+1+self.nz)
## Tgo to soil
self.dataTemp = np.ones(self.nz)*-self.F3F1
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+4*self.nz)
## Tgo to tout
self.dataTemp = np.ones(self.nz)*-self.F3F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+3*self.nz)
## Tgo to Tgi
self.dataTemp = np.ones(self.nz)*-2*self.F3F4
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
##### Tfi
## main diagonal Tfi
self.dataTemp = np.ones(self.nz)*(1+self.F5+self.F6F2)
self.dataTemp[0] = 1
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 2*self.nz)
## Tfi z-1
self.data = np.append(self.data,np.ones(self.nz-1)*-self.F5)
self.posX = np.append(self.posX,np.arange(0,self.nz-1,1)+1+ 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz-1,1)+ 2*self.nz)
## Tfi to Tgi
self.dataTemp = np.ones(self.nz)*-self.F6F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+ 2*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1))
##### Tfo
## main diagonal Tfo
self.dataTemp = np.ones(self.nz)*(1+self.F5+self.F6F2)
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 3*self.nz)
## Tfo z+1
self.data = np.append(self.data,np.ones(self.nz)*-self.F5)
self.posXTemp = np.arange(0,self.nz,1)+ 3*self.nz
self.posYTemp = np.arange(0,self.nz,1)+1+3*self.nz
self.posYTemp[self.nz-1] = self.posYTemp[self.nz-1]-(self.nz+1)
self.posX = np.append(self.posX,self.posXTemp)
self.posY = np.append(self.posY,self.posYTemp)
## Tfo to Tgo
self.dataTemp = np.ones(self.nz)*-self.F6F2
self.dataTemp[0] = 0
self.data = np.append(self.data,self.dataTemp)
self.posX = np.append(self.posX,np.arange(0,self.nz,1)+ 3*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1)+self.nz)
# main diagonal Soil 1
self.data = np.append(self.data,np.ones(self.nz))
self.posX = np.append(self.posX,np.arange(0,self.nz,1) + 4*self.nz)
self.posY = np.append(self.posY,np.arange(0,self.nz,1) + 4*self.nz)
self.K_NoFlow = csr_matrix((self.data, (self.posX, self.posY)), shape=(self.nz*self.na, self.nz*self.na),dtype=np.float)
self.K_sparse_NoFlow = sparse.csr_matrix(self.K_NoFlow)
######## Talt ########
######################
self.result = np.ones(self.nz*self.na)*BheData['Tundist'] # Initialbedingung
self.Tf_in = np.full(self.nz,BheData['Tundist']) # CellArray Tf_in
self.Tf_out = np.full(self.nz,BheData['Tundist']) # CellArray Tf_out
self.T_grout_in = np.full(self.nz,BheData['Tundist']) # CellArray Tf_grout
self.T_grout_out = np.full(self.nz,BheData['Tundist']) # CellArray Tf_grout
return True
def calcSondeFlow(self,tfinal,T_in):
for i in range(0,tfinal):
self.result[2*self.nz] = T_in # Vorgabe Randbedingung Tin
self.U =sparse.linalg.spsolve(self.K_sparse_Flow,self.result) # Lösen GlS
# Rausschreiben für Plot
self.T_grout_in = self.U[0:self.nz]
self.T_grout_out = self.U[self.nz:2*self.nz]
self.Tf_in = self.U[2*self.nz:3*self.nz]
self.Tf_out = self.U[3*self.nz:4*self.nz]
# reinschreiben nächster Zeitschritt
self.result[:] = np.copy(self.U[:])
return True
def calcSondeFlowQ(self,tfinal,T_in,Q):
if Q != self.Qold:
#if Q > self.Qold*1.05 or Q < self.Qold*0.95:
#if True:
self.Qold = Q
# Flow velocity
self.u = self.Qold/(2*np.pi*self.rpi**2)
# Flow parameters
self.Pr = self.dynviscF*self.capF/self.densF/self.lmF
self.Re = self.u*self.dpi/(self.dynviscF/self.densF)
self.Nu = Nusselt2U (self.Re,self.Pr,self.rpi,self.length)
self.Radv = 1/(self.Nu*self.lmF*np.pi)
self.Rfg = self.Radv + self.Rcona + self.Rconb
# Variables Flow
self.F2 = self.dz/self.Rfg
self.F5 = self.dt/self.dz*self.u
self.F6F2 = self.F6*self.F2
self.F3F2 = self.F3*self.F2
self.flow_vari[0] = (1 + self.F5 + self.F6F2)
self.flow_vari[1] = (- self.F5)
self.flow_vari[2] = (- self.F6F2)
self.flow_vari[3] = (1 + 2*self.F3F8 + self.F3F1 + self.F3F2 + 2*self.F3F4)
self.flow_vari[4] = (1 + self.F3F8 + self.F3F1 + self.F3F2 + 2*self.F3F4)
self.flow_vari[5] = (- self.F3F2)
# Write updated Variables to matrix
for i in range(0,6):
self.vari = self.vari_indices[i]
self.K_arr[tuple(self.vari)] = self.flow_vari[i]
self.K_sparse_Flow = sparse.csr_matrix(self.K_arr)
for i in range(0,tfinal):
self.result[2*self.nz] = T_in # Vorgabe Randbedingung Tin
self.U =sparse.linalg.spsolve(self.K_sparse_Flow,self.result) # Lösen GlS
# Rausschreiben für Plot
self.T_grout_in = self.U[0:self.nz]
self.T_grout_out = self.U[self.nz:2*self.nz]
self.Tf_in = self.U[2*self.nz:3*self.nz]
self.Tf_out = self.U[3*self.nz:4*self.nz]
# reinschreiben nächster Zeitschritt
self.result[:] = np.copy(self.U[:])
return True
def calcSondeNoFlow(self,tfinal):
for i in range(0,tfinal):
self.U = sparse.linalg.spsolve(self.K_sparse_NoFlow,self.result)
# Rausschreiben für Plot
self.T_grout_in = self.U[0:self.nz]
self.T_grout_out = self.U[self.nz:2*self.nz]
self.Tf_in = self.U[2*self.nz:3*self.nz]
self.Tf_out = self.U[3*self.nz:4*self.nz]
# reinschreiben nächster Zeitschritt
self.result[:] = np.copy(self.U[:])
return True
class BHE_2U_expl:
'''
explicit simulation model for 2U bhe
'''
def setTimestep(self,dt):
self.dt = dt
def setnz(self,nz):
self.nz = nz
def setSoilBC(self,Tsoil):
self.TsA[:] = Tsoil
def getGroutBC(self):
return (np.mean(self.T_grout_in)+np.mean(self.T_grout_out))/2
def getFluidOut(self):
return self.Tf_out[0]
def getFluidIn(self):
return self.Tf_in[0]
def initialize(self,BheData):
self.dz = BheData['length']/self.nz
self.dynviscF = BheData['dynviscF']
# Geometry
self.D = BheData['diamB']
self.ro = BheData['diamB']/2
self.rpi = BheData['odiamP']/2-BheData['thickP']
self.rpo = BheData['odiamP']/2
self.dpi = BheData['odiamP']-2*BheData['thickP']
self.dpo = BheData['odiamP']
# fluid properties
self.dynviscF = BheData['dynviscF']
self.densF = BheData['densF']
self.lmF = BheData['lmF']
self.length = BheData['length']
self.capF = BheData['capF']
# Flow velocity
self.u = BheData['Qf']/(2*np.pi*self.rpi**2)
self.Qold = BheData['Qf']
self.maxVel = self.u
# Flow parameters
self.Pr = self.dynviscF*BheData['capF']/BheData['densF']/BheData['lmF']
self.Re = self.u*self.dpi/(self.dynviscF/BheData['densF'])
self.Nu = Nusselt2U (self.Re,self.Pr,self.rpi,BheData['length'])
self.Radv = 1/(self.Nu*BheData['lmF']*np.pi)
# Resistances