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
|
/* Copyright 2019-2020 Yinjian Zhao
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include "FieldEnergy.H"
#include "WarpX.H"
#include "WarpXConst.H"
#include "AMReX_REAL.H"
#include "AMReX_ParticleReduce.H"
#include <iostream>
#include <cmath>
using namespace amrex;
// constructor
FieldEnergy::FieldEnergy (std::string rd_name)
: ReducedDiags{rd_name}
{
// RZ coordinate is not working
#if (defined WARPX_DIM_RZ)
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false,
"FieldEnergy reduced diagnostics does not work for RZ coordinate.");
#endif
// get WarpX class object
auto & warpx = WarpX::GetInstance();
// read number of levels
int nLevel = 0;
ParmParse pp("amr");
pp.query("max_level", nLevel);
nLevel += 1;
// resize data array
m_data.resize(3*nLevel,0.0);
if (ParallelDescriptor::IOProcessor())
{
if ( m_IsNotRestart )
{
// open file
std::ofstream ofs;
ofs.open(m_path + m_rd_name + "." + m_extension,
std::ofstream::out | std::ofstream::app);
// write header row
ofs << "#";
ofs << "[1]step";
ofs << m_sep;
ofs << "[2]time(s)";
for (int lev = 0; lev < nLevel; ++lev)
{
ofs << m_sep;
ofs << "[" + std::to_string(3+3*lev) + "]";
ofs << "total(J)lev"+std::to_string(lev);
ofs << m_sep;
ofs << "[" + std::to_string(4+3*lev) + "]";
ofs << "E(J)lev"+std::to_string(lev);
ofs << m_sep;
ofs << "[" + std::to_string(5+3*lev) + "]";
ofs << "B(J)lev"+std::to_string(lev);
}
ofs << std::endl;
// close file
ofs.close();
}
}
}
// end constructor
// function that computes field energy
void FieldEnergy::ComputeDiags (int step)
{
// Judge if the diags should be done
if ( (step+1) % m_freq != 0 ) { return; }
// get WarpX class object
auto & warpx = WarpX::GetInstance();
// get number of level
auto nLevel = warpx.finestLevel() + 1;
// loop over refinement levels
for (int lev = 0; lev < nLevel; ++lev)
{
// get MultiFab data at lev
const MultiFab & Ex = warpx.getEfield(lev,0);
const MultiFab & Ey = warpx.getEfield(lev,1);
const MultiFab & Ez = warpx.getEfield(lev,2);
const MultiFab & Bx = warpx.getBfield(lev,0);
const MultiFab & By = warpx.getBfield(lev,1);
const MultiFab & Bz = warpx.getBfield(lev,2);
// get cell size
Geometry const & geom = warpx.Geom(lev);
auto domain_box = geom.Domain();
#if (AMREX_SPACEDIM == 2)
auto dV = geom.CellSize(0) * geom.CellSize(1);
#elif (AMREX_SPACEDIM == 3)
auto dV = geom.CellSize(0) * geom.CellSize(1) * geom.CellSize(2);
#endif
// compute E squared
Real tmpx = Ex.norm2(0,geom.periodicity());
Real tmpy = Ey.norm2(0,geom.periodicity());
Real tmpz = Ez.norm2(0,geom.periodicity());
Real Es = tmpx*tmpx + tmpy*tmpy + tmpz*tmpz;
// compute B squared
tmpx = Bx.norm2(0,geom.periodicity());
tmpy = By.norm2(0,geom.periodicity());
tmpz = Bz.norm2(0,geom.periodicity());
Real Bs = tmpx*tmpx + tmpy*tmpy + tmpz*tmpz;
// save data
m_data[lev*3+1] = 0.5 * Es * PhysConst::ep0 * dV;
m_data[lev*3+2] = 0.5 * Bs / PhysConst::mu0 * dV;
m_data[lev*3+0] = m_data[lev*3+1] + m_data[lev*3+2];
}
// end loop over refinement levels
/* m_data now contains up-to-date values for:
* [total field energy at level 0,
* electric field energy at level 0,
* magnetic field energy at level 0,
* total field energy at level 1,
* electric field energy at level 1,
* magnetic field energy at level 1,
* ......] */
}
// end void FieldEnergy::ComputeDiags
|