-
Notifications
You must be signed in to change notification settings - Fork 9
/
sd1d.cxx
2283 lines (1870 loc) · 81.1 KB
/
sd1d.cxx
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
/*
SD1D: 1D simulation of plasma-neutral interactions
==================================================
Copyright B.Dudson, University of York, 2016-2018
email: [email protected]
This file is part of SD1D.
SD1D is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
SD1D is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with SD1D. If not, see <http://www.gnu.org/licenses/>.
Normalisations
--------------
Ne (density) normalised to Nnorm [m^-3]
T (temperature) normalised to Tnorm [eV]
B (magnetic field) normalised to Bnorm [eV]
t (time) normalised using ion cyclotron frequency Omega_ci [1/s]
Vi (velocity) normalised to sound speed Cs [m/s]
L (lengths) normalised to hybrid Larmor radius rho_s = Cs/Omega_ci [m]
*/
#include <mpi.h>
#include "revision.hxx"
#include <bout/constants.hxx>
#include <bout/physicsmodel.hxx>
#include <derivs.hxx>
#include <field_factory.hxx>
#include <bout/invert_pardiv.hxx>
#include <bout/snb.hxx>
#include <bout/fv_ops.hxx>
#include "div_ops.hxx"
#include "loadmetric.hxx"
#include "radiation.hxx"
// OpenADAS interface Atomicpp by T.Body
#include "atomicpp/ImpuritySpecies.hxx"
#include "atomicpp/Prad.hxx"
using bout::HeatFluxSNB;
// Helper for outputting variables
template<typename T>
void set_with_attrs(Options& option, T value, std::initializer_list<std::pair<std::string, Options::AttributeType>> attrs) {
option.force(value);
option.setAttributes(attrs);
}
class SD1D : public PhysicsModel {
protected:
int init(bool restarting) override {
Options &opt = Options::root()["sd1d"];
output.write("\nGit Version of SD1D: %s\n", sd1d::version::revision);
opt["revision"] = sd1d::version::revision;
opt["revision"].setConditionallyUsed();
OPTION(opt, cfl_info, false); // Calculate and print CFL information
// Normalisation
OPTION(opt, Tnorm, 100); // Reference temperature [eV]
OPTION(opt, Nnorm, 1e19); // Reference density [m^-3]
OPTION(opt, Bnorm, 1.0); // Reference magnetic field [T]
OPTION(opt, AA, 2.0); // Ion mass
// Model parameters
OPTION(opt, vwall, 1. / 3); // 1/3rd Franck-Condon energy at wall
OPTION(opt, frecycle, 1.0); // Recycling fraction 100%
OPTION(opt, fredistribute,
0.0); // Fraction of neutrals redistributed evenly along leg
OPTION(opt, density_sheath, 0); // Free boundary
OPTION(opt, pressure_sheath, 0); // Free boundary
OPTION(opt, gaspuff, 0.0); // Additional gas flux at target
OPTION(opt, include_dneut, true); // Include neutral gas diffusion?
if (opt.isSet("dneut")) {
// Scale neutral gas diffusion
dneut = FieldFactory::get()->create3D("dneut", &opt, mesh);
} else {
dneut = 1.0;
}
if (min(dneut) < 0.0) {
throw BoutException("dneut must be >= 0. Set include_dneut=false to disable\n");
}
OPTION(opt, nloss, 0.0); // Neutral gas loss rate
OPTION(opt, Eionize, 30); // Energy loss per ionisation (30eV)
OPTION(opt, sheath_gamma, 6.5); // Sheath heat transmission
OPTION(opt, neutral_gamma, 0.25); // Neutral heat transmission
// Plasma anomalous transport
OPTION(opt, anomalous_D, -1);
OPTION(opt, anomalous_chi, -1);
if (sheath_gamma < 6)
throw BoutException("sheath_gamma < 6 not consistent");
OPTION(opt, tn_floor, 3.5); // Minimum neutral gas temperature [eV]
OPTION(opt, atomic, true);
OPTION(opt, neutral_f_pn, true);
OPTION(opt, hyper, -1); // Numerical hyper-diffusion
OPTION(opt, ADpar, -1); // Added Dissipation scheme
OPTION(opt, viscos, -1); // Parallel viscosity
OPTION(opt, ion_viscosity, false); // Braginskii parallel viscosity
OPTION(opt, heat_conduction, true); // Spitzer-Hahm heat conduction
kappa_limit_alpha = opt["kappa_limit_alpha"]
.doc("Flux limiter. Turned off if < 0 (default)")
.withDefault(-1.0);
snb_model = opt["snb_model"]
.doc("Use SNB non-local heat flux model")
.withDefault<bool>(false);
if (snb_model) {
// Create a solver to calculate the SNB heat flux
snb = new HeatFluxSNB();
}
OPTION(opt, charge_exchange, true);
OPTION(opt, charge_exchange_escape, false);
OPTION(opt, charge_exchange_return_fE, 1.0);
OPTION(opt, recombination, true);
OPTION(opt, ionisation, true);
OPTION(opt, elastic_scattering,
false); // Include ion-neutral elastic scattering?
OPTION(opt, excitation, false); // Include electron impact excitation?
OPTION(opt, gamma_sound, 5. / 3); // Ratio of specific heats
bndry_flux_fix =
opt["bndry_flux_fix"]
.doc("Calculate boundary fluxes using simple mid-point (recommended)")
.withDefault<bool>(true);
// Field factory for generating fields from strings
FieldFactory ffact(mesh);
// Calculate normalisation factors
Cs0 = sqrt(SI::qe * Tnorm / (AA * SI::Mp)); // Reference sound speed [m/s]
Omega_ci = SI::qe * Bnorm / (AA * SI::Mp); // Ion cyclotron frequency [1/s]
rho_s0 = Cs0 / Omega_ci; // Length scale [m]
mi_me = AA * SI::Mp / SI::Me;
BoutReal Coulomb = 6.6 - 0.5 * log(Nnorm * 1e-20) + 1.5 * log(Tnorm);
tau_e0 = 1. / (2.91e-6 * (Nnorm / 1e6) * Coulomb * pow(Tnorm, -3. / 2));
OPTION(opt, volume_source, true);
OPTION(opt, time_dependent_source, false);
if (volume_source) {
// Volume sources of particles and energy
NeSource = Options::root()["Ne"]["source"]
.doc("Source of electron density. In SI units of particles/m^3/s").as<Field2D>();
if (time_dependent_source){
// time depentend sources Ne and P sources for elms
NeElm = FieldFactory::get()->create3D("Ne:NeElm");
ElmPrefactorGenerator = FieldFactory::get()->parse("Ne:ElmPrefactor");
PElm = FieldFactory::get()->create3D("P:PElm");
}
PeSource = Options::root()["P"]["source"]
.doc("Source of pressure in SI units of Pascals/s. Multiply by 3/2 to get W/m3/s").as<Field2D>();
// If the mesh file contains a source_weight variable, scale sources
Field2D source_weight;
if (mesh->get(source_weight, "source_weight") == 0) {
// Does have the source function in the input
output_info.write("Multiplying density and pressure sources by source_weight from mesh\n");
NeSource *= source_weight;
PeSource *= source_weight;
}
// Normalise sources
NeSource /= Nnorm * Omega_ci;
PeSource /= SI::qe * Nnorm * Tnorm * Omega_ci;
} else {
// Point sources, fixing density and specifying energy flux
Options *optpe = Options::getRoot()->getSection("P");
OPTION(optpe, powerflux, 2e7); // Power flux in W/m^2
powerflux /=
rho_s0 * SI::qe * Tnorm * Nnorm * Omega_ci; // Normalised energy flux
}
/////////////////////////
// Density controller
OPTION(opt, density_upstream, -1); // Fix upstream density? [m^-3]
if (density_upstream > 0.0) {
// Fixing density
density_upstream /= Nnorm;
// Controller
OPTION(opt, density_controller_p, 1e-2);
OPTION(opt, density_controller_i, 1e-3);
OPTION(opt, density_integral_positive, false);
OPTION(opt, density_source_positive, true);
density_error_lasttime = -1.0; // Signal no value
if (!restarting) {
density_error_integral = 0.0;
if (volume_source) {
// Set density_error_integral so that
// the input source is used
density_error_integral = 1. / density_controller_i;
}
}
}
if (volume_source) {
NeSource0 = NeSource; // Save initial value
}
Options::getRoot()->getSection("NVn")->get("evolve", evolve_nvn, true);
Options::getRoot()->getSection("Pn")->get("evolve", evolve_pn, true);
nloss /= Omega_ci;
// Specify variables to evolve
solver->add(Ne, "Ne");
solver->add(NVi, "NVi");
solver->add(P, "P");
if (atomic) {
solver->add(Nn, "Nn");
// Get the neutral density source
NnSource = Options::root()["Nn"]["source"]
.doc("Neutral atom source. SI units of particles/m^3/s").withDefault(Field2D(0.0))
/ (Nnorm * Omega_ci);
if (evolve_nvn) {
solver->add(NVn, "NVn");
}
if (evolve_pn) {
solver->add(Pn, "Pn");
// Get the neutral pressure source
PnSource = Options::root()["Pn"]["source"]
.doc("Neutral atom pressure source. SI units of Pa/s").withDefault(Field2D(0.0))
/ (SI::qe * Nnorm * Tnorm * Omega_ci);
}
}
// Load the metric tensor
LoadMetric(rho_s0, Bnorm);
if ( opt.isSet("area") ) {
// Area set in the input file. Overwrite any Jacobian from the mesh
mesh->getCoordinates()->J = ffact.create2D(opt["area"].as<std::string>(),
Options::getRoot());
}
dy4 = SQ(SQ(mesh->getCoordinates()->dy));
//////////////////////////////////////////////////
// Impurities
OPTION(opt, fimp, 0.0); // Fixed impurity fraction
OPTION(opt, impurity_adas, false);
if (impurity_adas) {
// Use OpenADAS data through Atomicpp
// Find out which species to model
string impurity_species;
OPTION(opt, impurity_species, "c");
impurity = new ImpuritySpecies(impurity_species);
} else {
// Use carbon radiation for the impurity
rad = new HutchinsonCarbonRadiation();
}
diagnose = opt["diagnose"].doc("Output extra variables").withDefault<bool>(true);
output_ddt =
opt["output_ddt"].doc("Save time derivatives?").withDefault<bool>(false);
kappa_epar = 0.0;
Srec = 0.0;
Siz = 0.0;
S = 0.0;
Frec = 0.0;
Fiz = 0.0;
Fcx = 0.0;
Fel = 0.0;
F = 0.0;
Rrec = 0.0;
Riz = 0.0;
Rzrad = 0.0;
Rex = 0.0;
R = 0.0;
Erec = 0.0;
Eiz = 0.0;
Ecx = 0.0;
Eel = 0.0;
E = 0.0;
Dcx = 0.0;
Dcx_T = 0.0;
flux_ion = 0.0;
// Neutral gas diffusion and heat conduction
Dn = 0.0;
kappa_n = 0.0;
// Anomalous transport
if (anomalous_D > 0.0) {
// Normalise
anomalous_D /= rho_s0 * rho_s0 * Omega_ci; // m^2/s
output.write("\tnormalised anomalous D_perp = %e\n", anomalous_D);
}
if (anomalous_chi > 0.0) {
// Normalise
anomalous_chi /= rho_s0 * rho_s0 * Omega_ci; // m^2/s
output.write("\tnormalised anomalous chi_perp = %e\n", anomalous_chi);
}
// Calculate neutral gas redistribution weights over the domain
string redist_string;
opt.get("redist_weight", redist_string, "1.0");
redist_weight = ffact.create2D(redist_string, &opt);
BoutReal localweight = 0.0;
Coordinates *coord = mesh->getCoordinates();
for (int j = mesh->ystart; j <= mesh->yend; j++) {
localweight += redist_weight(mesh->xstart, j) *
coord->J(mesh->xstart, j) * coord->dy(mesh->xstart, j);
}
MPI_Comm ycomm = mesh->getYcomm(mesh->xstart); // MPI communicator
// Calculate total weight by summing over all processors
BoutReal totalweight;
MPI_Allreduce(&localweight, &totalweight, 1, MPI_DOUBLE, MPI_SUM, ycomm);
// Normalise redist_weight so sum over domain:
//
// sum ( redist_weight * J * dy ) = 1
//
redist_weight /= totalweight;
setPrecon((preconfunc)&SD1D::precon);
//////////////////////////////////////////
// Split operator (IMEX) schemes
// Use combination of explicit and implicit methods
//
// Boolean flags rhs_explicit and rhs_implicit
// turn on explicit and implicit terms
bool split_operator;
OPTION(opt, split_operator, false);
if (!split_operator) {
// Turn on all terms in rhs
rhs_explicit = rhs_implicit = true;
update_coefficients = true;
}
setSplitOperator(split_operator);
return 0;
}
/*!
* This function calculates the time derivatives
* of all evolving quantities
*
*/
int rhs(BoutReal time) override {
Coordinates *coord = mesh->getCoordinates();
mesh->communicate(Ne, NVi, P);
// Floor small values
P = floor(P, 1e-10);
Ne = floor(Ne, 1e-10);
Field3D Nelim = floor(Ne, 1e-5);
Vi = NVi / Ne;
Field3D Te = 0.5 * P / Ne; // Assuming Te = Ti
for (auto &i : Te.getRegion("RGN_NOBNDRY")) {
if (Te[i] > 10.)
Te[i] = 10.;
}
Field3D Nnlim;
Field3D Tn;
if (atomic) {
// Includes atomic processes, neutral gas
mesh->communicate(Nn);
if (evolve_nvn) {
mesh->communicate(NVn);
}
if (evolve_pn) {
mesh->communicate(Pn);
}
Nn = floor(Nn, 1e-10);
Nnlim = floor(Nn, 1e-5);
if (evolve_nvn) {
Vn = NVn / Nnlim;
} else {
Vn = -vwall * sqrt(3.5 / Tnorm);
NVn = Nn * Vn;
}
if (evolve_pn) {
Tn = Pn / Nnlim;
// Tn = floor(Tn, 0.025/Tnorm); // Minimum tn_floor
Tn = floor(Tn, 1e-12);
} else {
Tn = Te; // Strong CX coupling
Pn = Tn * floor(Nn, 0.0);
Tn = floor(Tn, tn_floor / Tnorm); // Minimum of tn_floor
}
}
if (update_coefficients) {
// Update diffusion coefficients
TRACE("Update coefficients");
tau_e = Omega_ci * tau_e0 * pow(Te, 1.5) / Ne;
if (heat_conduction) {
kappa_epar = 3.2 * mi_me * 0.5 * P * tau_e;
if (kappa_limit_alpha > 0.0) {
/*
* Flux limiter, as used in SOLPS.
*
* Calculate the heat flux from Spitzer-Harm and flux limit
*
* Typical value of alpha ~ 0.2 for electrons
*
* R.Schneider et al. Contrib. Plasma Phys. 46, No. 1-2, 3 – 191 (2006)
* DOI 10.1002/ctpp.200610001
*/
// Spitzer-Harm heat flux
Te.applyBoundary("neumann"); // Note: We haven't yet applied boundaries
Field3D q_SH = kappa_epar * Grad_par(Te);
Field3D q_fl = kappa_limit_alpha * Nelim * Te * sqrt(mi_me * Te);
kappa_epar = kappa_epar / (1. + abs(q_SH / q_fl));
// Values of kappa on cell boundaries are needed for fluxes
mesh->communicate(kappa_epar);
}
kappa_epar.applyBoundary("neumann");
}
if (atomic) {
// Neutral diffusion rate
for (int i = 0; i < mesh->LocalNx; i++)
for (int j = 0; j < mesh->LocalNy; j++)
for (int k = 0; k < mesh->LocalNz; k++) {
// Charge exchange frequency, normalised to ion cyclotron
// frequency
BoutReal sigma_cx = Nelim(i, j, k) * Nnorm *
hydrogen.chargeExchange(Te(i, j, k) * Tnorm) /
Omega_ci;
// Ionisation frequency
BoutReal sigma_iz = Nelim(i, j, k) * Nnorm *
hydrogen.ionisation(Te(i, j, k) * Tnorm) /
Omega_ci;
// Neutral thermal velocity
BoutReal tn = Tn(i, j, k);
if (tn < tn_floor / Tnorm)
tn = tn_floor / Tnorm;
BoutReal vth_n = sqrt(tn); // Normalised to Cs0
// Neutral-neutral mean free path
BoutReal Lmax = 0.1; // meters
BoutReal a0 = PI * SQ(5.29e-11);
BoutReal lambda_nn = 1. / (Nnorm * Nnlim(i, j, k) * a0); // meters
if (lambda_nn > Lmax) {
// Limit maximum mean free path
lambda_nn = Lmax;
}
lambda_nn /= rho_s0; // Normalised length to rho_s0
// Neutral-Neutral collision rate
BoutReal sigma_nn = vth_n / lambda_nn;
// Total neutral collision frequency
BoutReal sigma = sigma_cx + sigma_iz + sigma_nn;
// Neutral gas diffusion
if (include_dneut) {
Dn(i, j, k) = dneut(i, j, k) * SQ(vth_n) / sigma;
// Neutral gas heat conduction
kappa_n(i, j, k) = dneut(i, j, k) * Nnlim(i, j, k) * SQ(vth_n) / sigma;
}
}
kappa_n.applyBoundary("Neumann");
Dn.applyBoundary("dirichlet_o2");
mesh->communicate(kappa_n, Dn);
}
}
// Set sheath boundary condition on flow
TRACE("Sheath");
ddt(P) = 0.0; // Need to set heat flux
if (evolve_pn) {
ddt(Pn) = 0.0;
}
for (RangeIterator r = mesh->iterateBndryUpperY(); !r.isDone(); r++) {
int jz = 0;
// Outward flow velocity to >= Cs
BoutReal Vout =
sqrt(2.0 * Te(r.ind, mesh->yend, jz)); // Sound speed outwards
if (Vi(r.ind, mesh->yend, jz) > Vout)
Vout = Vi(r.ind, mesh->yend,
jz); // If plasma is faster, go to plasma velocity
BoutReal Nout;
switch (density_sheath) {
case 0: {
// Free boundary on density (constant gradient)
Nout = 0.5 *
(3. * Ne(r.ind, mesh->yend, jz) - Ne(r.ind, mesh->yend - 1, jz));
break;
}
case 1: {
// Zero gradient
Nout = Ne(r.ind, mesh->yend, jz);
break;
}
case 2: {
// Zero gradient particle flux N*Vi* J*dx*dz
// Since Vi increases into the sheath, density should drop
Nout =
Ne(r.ind, mesh->yend, jz) * coord->J(r.ind, mesh->yend) *
Vi(r.ind, mesh->yend, jz) /
(0.5 *
(coord->J(r.ind, mesh->yend) + coord->J(r.ind, mesh->yend + 1)) *
Vout);
break;
}
default:
throw BoutException("Unrecognised density_sheath option");
}
if (Nout < 0.0)
Nout = 0.0; // Prevent Nout from going negative -> Flux is always to the
// wall
// Flux of particles is Ne*Vout
BoutReal flux = Nout * Vout;
BoutReal Pout;
switch (pressure_sheath) {
case 0: {
// Free boundary (constant gradient)
Pout = 0.5 *
(3. * P(r.ind, mesh->yend, jz) - P(r.ind, mesh->yend - 1, jz));
break;
}
case 1: {
// Zero gradient
Pout = P(r.ind, mesh->yend, jz);
break;
}
case 2: {
// Use energy flux conservation to set pressure
// (5/2)Pv + (1/2)nv^3 = const
//
Pout =
((5. * P(r.ind, mesh->yend, jz) * Vi(r.ind, mesh->yend, jz) +
Ne(r.ind, mesh->yend, jz) * pow(Vi(r.ind, mesh->yend, jz), 3)) /
Vout -
Nout * Vout * Vout) /
5.;
break;
}
default:
throw BoutException("Unrecognised pressure_sheath option");
}
if (Pout < 0.0)
Pout = 0.0;
if (rhs_explicit) {
// Additional cooling
BoutReal q = (sheath_gamma - 6) * Te(r.ind, mesh->yend, jz) * flux;
// Multiply by cell area to get power
BoutReal heatflux =
q *
(coord->J(r.ind, mesh->yend) + coord->J(r.ind, mesh->yend + 1)) /
(sqrt(coord->g_22(r.ind, mesh->yend)) +
sqrt(coord->g_22(r.ind, mesh->yend + 1)));
// Divide by volume of cell, and 2/3 to get pressure
ddt(P)(r.ind, mesh->yend, jz) -=
(2. / 3) * heatflux /
(coord->dy(r.ind, mesh->yend) * coord->J(r.ind, mesh->yend));
}
// Set boundary half-way between cells
for (int jy = mesh->yend + 1; jy < mesh->LocalNy; jy++) {
///// Plasma model
// Vi fixed value (Dirichlet)
Vi(r.ind, jy, jz) = 2. * Vout - Vi(r.ind, mesh->yend, jz);
// Ne set from flux (Dirichlet)
Ne(r.ind, jy, jz) = 2 * Nout - Ne(r.ind, mesh->yend, jz);
// NVi. This can be negative, so set this to the flux
// going out of the domain (zero gradient)
NVi(r.ind, jy, jz) = Nout * Vout;
// NVi(r.ind, jy, jz) = Ne(r.ind, jy, jz) * Vi(r.ind, jy, jz);
// NVi(r.ind, jy, jz) = 2.*Nout * Vout - NVi(r.ind, mesh->yend, jz);
// Te zero gradient (Neumann)
Te(r.ind, jy, jz) = Te(r.ind, mesh->yend, jz);
P(r.ind, jy, jz) = 2. * Pout - P(r.ind, mesh->yend, jz);
if (atomic) {
///// Neutral model
// Flux of neutral particles, momentum, and energy are set later
// Here the neutral velocity is set to no-flow conditions
// Vn fixed value (Dirichlet)
Vn(r.ind, jy, jz) = -Vn(r.ind, mesh->yend, jz);
// Nn free boundary (constant gradient)
Nn(r.ind, jy, jz) =
2. * Nn(r.ind, mesh->yend, jz) - Nn(r.ind, mesh->yend - 1, jz);
if (evolve_pn) {
// Tn fixed value (Dirichlet)
// Tn(r.ind, jy, jz) = 3.5/Tnorm - Tn(r.ind, mesh->yend, jz);
// Tn zero gradient. Heat flux set by gamma
Tn(r.ind, jy, jz) = Tn(r.ind, mesh->yend, jz);
if (rhs_explicit && (neutral_gamma > 0.0)) {
// Density at the target
BoutReal Nnout = 0.5 * (Nn(r.ind, mesh->yend, jz) +
Nn(r.ind, mesh->yend + 1, jz));
// gamma * n * T * cs
BoutReal q = neutral_gamma * Nnout * Tn(r.ind, jy, jz) *
sqrt(Tn(r.ind, jy, jz));
// Multiply by cell area to get power
BoutReal heatflux = q *
(coord->J(r.ind, mesh->yend) +
coord->J(r.ind, mesh->yend + 1)) /
(sqrt(coord->g_22(r.ind, mesh->yend)) +
sqrt(coord->g_22(r.ind, mesh->yend + 1)));
// Divide by volume of cell, and 2/3 to get pressure
ddt(Pn)(r.ind, mesh->yend, jz) -=
(2. / 3) * heatflux /
(coord->dy(r.ind, mesh->yend) * coord->J(r.ind, mesh->yend));
}
} else {
Tn(r.ind, jy, jz) = Te(r.ind, jy, jz);
}
Pn(r.ind, jy, jz) = Nn(r.ind, jy, jz) * Tn(r.ind, jy, jz);
NVn(r.ind, jy, jz) = -NVn(r.ind, mesh->yend, jz);
}
}
}
for (RangeIterator r = mesh->iterateBndryLowerY(); !r.isDone(); r++) {
// No-flow boundary condition on left boundary
for (int jz = 0; jz < mesh->LocalNz; jz++) {
for (int jy = 0; jy < mesh->ystart; jy++) {
Te(r.ind, jy, jz) = Te(r.ind, mesh->ystart, jz);
Ne(r.ind, jy, jz) = Ne(r.ind, mesh->ystart, jz);
P(r.ind, jy, jz) = P(r.ind, mesh->ystart, jz);
Vi(r.ind, jy, jz) = -Vi(r.ind, mesh->ystart, jz);
NVi(r.ind, jy, jz) = -NVi(r.ind, mesh->ystart, jz);
if (atomic) {
Vn(r.ind, jy, jz) = -Vn(r.ind, mesh->ystart, jz);
Nn(r.ind, jy, jz) = Nn(r.ind, jy, jz);
Pn(r.ind, jy, jz) = Pn(r.ind, jy, jz);
Tn(r.ind, jy, jz) = Tn(r.ind, jy, jz);
}
}
}
}
if ((density_upstream > 0.0) && rhs_explicit) {
///////////////////////////////////////////////
// Set velocity on left boundary to set density
//
// This calculates a source needed in the first grid cell, to relax
// towards the desired density value.
//
TRACE("Density upstream");
BoutReal source;
for (RangeIterator r = mesh->iterateBndryLowerY(); !r.isDone(); r++) {
int jz = 0;
// Density source, so dn/dt = source
BoutReal error = density_upstream - Ne(r.ind, mesh->ystart, jz);
ASSERT2(std::isfinite(error));
ASSERT2(std::isfinite(density_error_integral));
// PI controller, using crude integral of the error
if (density_error_lasttime < 0.0) {
// First time
density_error_lasttime = time;
density_error_last = error;
}
// Integrate using Trapezium rule
if (time > density_error_lasttime) { // Since time can decrease
density_error_integral += (time - density_error_lasttime) * 0.5 *
(error + density_error_last);
}
if ((density_error_integral < 0.0) && density_integral_positive) {
// Limit density_error_integral to be >= 0
density_error_integral = 0.0;
}
// Calculate source from combination of error and integral
source = density_controller_p * error +
density_controller_i * density_error_integral;
// output.write("\n Source: %e, %e : %e, %e -> %e\n", time, (time -
// density_error_lasttime), error, density_error_integral, source);
density_error_last = error;
density_error_lasttime = time;
if (!volume_source) {
// Convert source into a flow velocity
// through the boundary, based on a zero-gradient boundary on the
// density. This ensures that the mass and momentum inputs are
// consistent, but also carries energy through the boundary. This flux
// of energy is calculated, and subtracted from the pressure equation,
// so that the density boundary does not contribute to energy balance.
// Calculate needed input velocity
BoutReal Vin = source * sqrt(coord->g_22(r.ind, mesh->ystart)) *
coord->dy(r.ind, mesh->ystart) /
Ne(r.ind, mesh->ystart, jz);
// Limit at sound speed
BoutReal cs = sqrt(Te(r.ind, mesh->ystart, jz));
if (fabs(Vin) > cs) {
Vin *= cs / fabs(Vin); // + or - cs
}
Vi(r.ind, mesh->ystart - 1, jz) =
2. * Vin - Vi(r.ind, mesh->ystart, jz);
// Power flux is v * (5/2 P + 1/2 m n v^2 )
BoutReal inputflux = Vin * (2.5 * P(r.ind, mesh->ystart, jz) +
0.5 * Ne(r.ind, mesh->ystart, jz) * Vin *
Vin); // W/m^2 (normalised)
// Subtract input energy flux from P equation
// so no net power input
ddt(P)(r.ind, mesh->ystart, jz) -=
(2. / 3) * inputflux /
(coord->dy(r.ind, mesh->ystart) *
sqrt(coord->g_22(r.ind, mesh->ystart)));
}
}
if (volume_source) {
if ((source < 0.0) && density_source_positive) {
source = 0.0; // Don't remove particles
}
// Broadcast the value of source from processor 0
MPI_Bcast(&source, 1, MPI_DOUBLE, 0, BoutComm::get());
ASSERT2(std::isfinite(source));
// Scale NeSource
NeSource = source * NeSource0;
}
}
if (atomic && rhs_explicit) {
// Atomic physics
TRACE("Atomic");
// Lower floor on Nn for atomic rates
Field3D Nnlim2 = floor(Nn, 0.0);
if (fimp > 0.0) {
// Impurity radiation
if (impurity_adas) {
Rzrad.allocate();
for (auto &i : Rzrad.getRegion("RGN_NOY")) {
// Calculate cell centre (C), left (L) and right (R) values
BoutReal Te_C = Te[i],
Te_L = 0.5 * (Te[i.ym()] + Te[i]),
Te_R = 0.5 * (Te[i] + Te[i.yp()]);
BoutReal Ne_C = Ne[i],
Ne_L = 0.5 * (Ne[i.ym()] + Ne[i]),
Ne_R = 0.5 * (Ne[i] + Ne[i.yp()]);
BoutReal Nn_C = Nnlim2[i],
Nn_L = 0.5 * (Nnlim2[i.ym()] + Nnlim2[i]),
Nn_R = 0.5 * (Nnlim2[i] + Nnlim2[i.yp()]);
BoutReal Rz_L = computeRadiatedPower(*impurity,
Te_L * Tnorm, // electron temperature [eV]
Ne_L * Nnorm, // electron density [m^-3]
fimp * Ne_L * Nnorm, // impurity density [m^-3]
Nn_L * Nnorm); // Neutral density [m^-3]
BoutReal Rz_C = computeRadiatedPower(*impurity,
Te_C * Tnorm, // electron temperature [eV]
Ne_C * Nnorm, // electron density [m^-3]
fimp * Ne_C * Nnorm, // impurity density [m^-3]
Nn_C * Nnorm); // Neutral density [m^-3]
BoutReal Rz_R = computeRadiatedPower(*impurity,
Te_R * Tnorm, // electron temperature [eV]
Ne_R * Nnorm, // electron density [m^-3]
fimp * Ne_R * Nnorm, // impurity density [m^-3]
Nn_R * Nnorm); // Neutral density [m^-3]
// Jacobian (Cross-sectional area)
BoutReal J_C = coord->J[i],
J_L = 0.5 * (coord->J[i.ym()] + coord->J[i]),
J_R = 0.5 * (coord->J[i] + coord->J[i.yp()]);
// Simpson's rule, calculate average over cell
Rzrad[i] = (J_L * Rz_L +
4. * J_C * Rz_C +
J_R * Rz_R) / (6. * J_C);
}
} else {
Rzrad = rad->power(Te * Tnorm, Ne * Nnorm,
Ne * (Nnorm * fimp)); // J / m^3 / s
}
Rzrad /= SI::qe * Tnorm * Nnorm * Omega_ci; // Normalise
} // else Rzrad = 0.0 set in init()
E = 0.0; // Energy transfer to neutrals
for (int i = 0; i < mesh->LocalNx; i++)
for (int j = mesh->ystart; j <= mesh->yend; j++)
for (int k = 0; k < mesh->LocalNz; k++) {
// Integrate rates over each cell using Simpson's rule
// Calculate cell centre (C), left (L) and right (R) values
BoutReal Te_C = Te(i, j, k),
Te_L = 0.5 * (Te(i, j - 1, k) + Te(i, j, k)),
Te_R = 0.5 * (Te(i, j, k) + Te(i, j + 1, k));
BoutReal Ne_C = Ne(i, j, k),
Ne_L = 0.5 * (Ne(i, j - 1, k) + Ne(i, j, k)),
Ne_R = 0.5 * (Ne(i, j, k) + Ne(i, j + 1, k));
BoutReal Vi_C = Vi(i, j, k),
Vi_L = 0.5 * (Vi(i, j - 1, k) + Vi(i, j, k)),
Vi_R = 0.5 * (Vi(i, j, k) + Vi(i, j + 1, k));
BoutReal Tn_C = Tn(i, j, k),
Tn_L = 0.5 * (Tn(i, j - 1, k) + Tn(i, j, k)),
Tn_R = 0.5 * (Tn(i, j, k) + Tn(i, j + 1, k));
BoutReal Nn_C = Nnlim2(i, j, k),
Nn_L = 0.5 * (Nnlim2(i, j - 1, k) + Nnlim2(i, j, k)),
Nn_R = 0.5 * (Nnlim2(i, j, k) + Nnlim2(i, j + 1, k));
BoutReal Vn_C = Vn(i, j, k),
Vn_L = 0.5 * (Vn(i, j - 1, k) + Vn(i, j, k)),
Vn_R = 0.5 * (Vn(i, j, k) + Vn(i, j + 1, k));
// Jacobian (Cross-sectional area)
BoutReal J_C = coord->J(i, j),
J_L = 0.5 * (coord->J(i, j - 1) + coord->J(i, j)),
J_R = 0.5 * (coord->J(i, j) + coord->J(i, j + 1));
///////////////////////////////////////
// Charge exchange
if (charge_exchange) {
BoutReal R_cx_L = Ne_L * Nn_L *
hydrogen.chargeExchange(Te_L * Tnorm) *
(Nnorm / Omega_ci);
BoutReal R_cx_C = Ne_C * Nn_C *
hydrogen.chargeExchange(Te_C * Tnorm) *
(Nnorm / Omega_ci);
BoutReal R_cx_R = Ne_R * Nn_R *
hydrogen.chargeExchange(Te_R * Tnorm) *
(Nnorm / Omega_ci);
// Ecx is energy transferred to neutrals
Ecx(i, j, k) = (3. / 2) *
(J_L * (Te_L - Tn_L) * R_cx_L +
4. * J_C * (Te_C - Tn_C) * R_cx_C +
J_R * (Te_R - Tn_R) * R_cx_R) /
(6. * J_C);
// Fcx is friction between plasma and neutrals
Fcx(i, j, k) = (J_L * (Vi_L - Vn_L) * R_cx_L +
4. * J_C * (Vi_C - Vn_C) * R_cx_C +
J_R * (Vi_R - Vn_R) * R_cx_R) /
(6. * J_C);
// Dcx is a redistribution of fast neutrals due to charge exchange
// Acts as a sink of plasma density
Dcx(i, j, k) = (J_L * R_cx_L + 4. * J_C * R_cx_C + J_R * R_cx_R) /
(6. * J_C);
// Energy lost from the plasma
// This gives the temperature of the CX neutrals when
// divided by Dcx
Dcx_T(i, j, k) = (J_L * Te_L * R_cx_L + 4. * J_C * Te_C * R_cx_C +
J_R * Te_R * R_cx_R) /
(6. * J_C);
}
///////////////////////////////////////
// Recombination
if (recombination) {
BoutReal R_rc_L =
hydrogen.recombination(Ne_L * Nnorm, Te_L * Tnorm) *
SQ(Ne_L) * Nnorm / Omega_ci;
BoutReal R_rc_C =
hydrogen.recombination(Ne_C * Nnorm, Te_C * Tnorm) *
SQ(Ne_C) * Nnorm / Omega_ci;
BoutReal R_rc_R =
hydrogen.recombination(Ne_R * Nnorm, Te_R * Tnorm) *
SQ(Ne_R) * Nnorm / Omega_ci;
// Rrec is radiated energy, Erec is energy transferred to neutrals
// Factor of 1.09 so that recombination becomes an energy source
// at 5.25eV
Rrec(i, j, k) =
(J_L * (1.09 * Te_L - 13.6 / Tnorm) * R_rc_L +
4. * J_C * (1.09 * Te_C - 13.6 / Tnorm) * R_rc_C +
J_R * (1.09 * Te_R - 13.6 / Tnorm) * R_rc_R) /
(6. * J_C);
Erec(i, j, k) = (3. / 2) *
(J_L * Te_L * R_rc_L + 4. * J_C * Te_C * R_rc_C +
J_R * Te_R * R_rc_R) /
(6. * J_C);
Frec(i, j, k) = (J_L * Vi_L * R_rc_L + 4. * J_C * Vi_C * R_rc_C +
J_R * Vi_R * R_rc_R) /
(6. * J_C);
Srec(i, j, k) =
(J_L * R_rc_L + 4. * J_C * R_rc_C + J_R * R_rc_R) /
(6. * J_C);
}
///////////////////////////////////////
// Ionisation
if (ionisation) {
BoutReal R_iz_L = Ne_L * Nn_L *
hydrogen.ionisation(Te_L * Tnorm) * Nnorm /
Omega_ci;
BoutReal R_iz_C = Ne_C * Nn_C *
hydrogen.ionisation(Te_C * Tnorm) * Nnorm /