Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Revamped add fields::dft_fields_norm2 #1740

Merged
merged 12 commits into from
Sep 17, 2021
77 changes: 5 additions & 72 deletions python/adjoint/optimization_problem.py
Original file line number Diff line number Diff line change
Expand Up @@ -69,7 +69,7 @@ def __init__(
fcen=None,
df=None,
nf=None,
decay_dt=50,
decay_dt=25,
decay_fields=[mp.Ez],
decay_by=1e-6,
decimation_factor=1,
Expand Down Expand Up @@ -237,12 +237,8 @@ def forward_run(self):
# Forward run
self.sim.run(until_after_sources=stop_when_dft_decayed(
self.sim,
self.design_region_monitors,
self.decay_dt,
self.decay_fields,
self.fcen_idx,
self.decay_by,
True,
self.minimum_run_time,
self.maximum_run_time,
))
Expand Down Expand Up @@ -312,12 +308,8 @@ def adjoint_run(self):
# Adjoint run
self.sim.run(until_after_sources=stop_when_dft_decayed(
self.sim,
self.design_region_monitors,
self.decay_dt,
self.decay_fields,
self.fcen_idx,
self.decay_by,
True,
self.minimum_run_time,
self.maximum_run_time,
))
Expand Down Expand Up @@ -431,12 +423,8 @@ def calculate_fd_gradient(

self.sim.run(until_after_sources=stop_when_dft_decayed(
self.sim,
self.forward_monitors,
self.decay_dt,
self.decay_fields,
self.fcen_idx,
self.decay_by,
True,
self.minimum_run_time,
self.maximum_run_time,
))
Expand Down Expand Up @@ -466,12 +454,8 @@ def calculate_fd_gradient(
# add monitor used to track dft convergence
self.sim.run(until_after_sources=stop_when_dft_decayed(
self.sim,
self.forward_monitors,
self.decay_dt,
self.decay_fields,
self.fcen_idx,
self.decay_by,
True,
self.minimum_run_time,
self.maximum_run_time,
))
Expand Down Expand Up @@ -533,32 +517,16 @@ def plot2D(self, init_opt=False, **kwargs):

def stop_when_dft_decayed(
simob,
mon,
dt,
c,
fcen_idx,
decay_by,
yee_grid=False,
minimum_run_time=0,
maximum_run_time=None,
):
'''Step function that monitors the relative change in DFT fields for a list of monitors.

mon ............. a list of monitors
c ............... a list of components to monitor

'''
# get monitor dft output array dimensions
dims = []
for m in mon:
ci_dims = []
for ci in c:
comp = ci if yee_grid else mp.Dielectric
ci_dims += [simob.get_array_slice_dimensions(comp, vol=m.where)[0]]
dims.append(ci_dims)

# Record data in closure so that we can persitently edit
closure = {'previous_fields': [[None for di in d] for d in dims], 't0': 0}
closure = {'previous_fields': 1, 't0': 0}

def _stop(sim):
if sim.round_time() <= dt + closure['t0']:
Expand All @@ -567,44 +535,9 @@ def _stop(sim):
return True
else:
previous_fields = closure['previous_fields']

# Pull all the relevant frequency and spatial dft points
relative_change = []
current_fields = [[0 for di in d] for d in dims]
for mi, m in enumerate(mon):
for ic, cc in enumerate(c):
if isinstance(m, mp.DftFlux):
current_fields[mi][ic] = mp.get_fluxes(m)[fcen_idx]
elif isinstance(m, mp.DftFields):
current_fields[mi][ic] = atleast_3d(
sim.get_dft_array(m, cc, fcen_idx))
elif isinstance(m, mp.simulation.DftNear2Far):
current_fields[mi][ic] = atleast_3d(
sim.get_dft_array(m, cc, fcen_idx))
else:
raise TypeError(
"Monitor of type {} not supported".format(type(m)))

if previous_fields[mi][ic] is not None:
if np.sum(
np.abs(previous_fields[mi][ic] -
current_fields[mi][ic])) == 0:
relative_change.append(0)
elif np.all(np.abs(previous_fields[mi][ic])):
relative_change_raw = np.abs(
previous_fields[mi][ic] -
current_fields[mi][ic]) / np.abs(
previous_fields[mi][ic])
relative_change.append(
np.mean(relative_change_raw.flatten())
) # average across space and frequency
else:
relative_change.append(1)
else:
relative_change.append(1)

relative_change = np.mean(
relative_change) # average across monitors
current_fields = simob.fields.dft_fields_norm()
relative_change = np.abs(previous_fields-current_fields) / previous_fields
stevengj marked this conversation as resolved.
Show resolved Hide resolved

closure['previous_fields'] = current_fields
closure['t0'] = sim.round_time()

Expand Down
58 changes: 58 additions & 0 deletions src/dft.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -265,6 +265,64 @@ void dft_chunk::update_dft(double time) {
}
}

/* Return the L2 norm of all of the fields used to update DFTs. This is useful
to check whether the simulation is finished (whether all relevant fields have decayed).
(Collective operation.) */
double fields::dft_fields_norm() {
am_now_working_on(Other);
double sum = 0.0;
for (int i = 0; i < num_chunks; i++)
if (chunks[i]->is_mine()) sum += chunks[i]->dft_fields_norm2();
finished_working();
return std::sqrt(sum_to_all(sum));
}

double fields_chunk::dft_fields_norm2() const {
double sum = 0.0;
for (dft_chunk *cur = dft_chunks; cur; cur = cur->next_in_chunk)
sum += cur->dft_fields_norm2();
return sum;
}

static double sqr(std::complex<realnum> x) { return (x*std::conj(x)).real(); }

double dft_chunk::dft_fields_norm2() const {
if (!fc->f[c][0]) return 0.0;
int numcmp = fc->f[c][1] ? 2 : 1;
double sum = 0.0;
size_t idx_dft = 0;
const int Nomega = omega.size();
LOOP_OVER_IVECS(fc->gv, is, ie, idx) {
for (int i = 0; i < Nomega; ++i)
sum += sqr(dft[Nomega * idx_dft + i]);
idx_dft++;
}
return sum;
}

// return the minimum abs(freq) over all DFT chunks
stevengj marked this conversation as resolved.
Show resolved Hide resolved
double fields::dft_minfreq() const {
double minfreq = meep::infinity;
for (int i = 0; i < num_chunks; i++)
if (chunks[i]->is_mine())
minfreq = std::min(minfreq, chunks[i]->dft_minfreq());
return -max_to_all(-minfreq); // == min_to_all(minfreq)
}

double fields_chunk::dft_minfreq() const {
double minomega = meep::infinity;
for (dft_chunk *cur = dft_chunks; cur; cur = cur->next_in_chunk)
minomega = std::min(minomega, cur->minomega());
return minomega / (2*meep::pi);
}

double dft_chunk::minomega() const {
double minomega = meep::infinity;
for (const auto& o : omega)
minomega = std::min(minomega, std::abs(o));
return minomega;
}

void dft_chunk::scale_dft(complex<double> scale) {
for (size_t i = 0; i < N * omega.size(); ++i)
dft[i] *= scale;
Expand Down
6 changes: 6 additions & 0 deletions src/meep.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1110,6 +1110,8 @@ class dft_chunk {
~dft_chunk();

void update_dft(double time);
double dft_fields_norm2() const;
double minomega() const;

void scale_dft(std::complex<double> scale);

Expand Down Expand Up @@ -1561,6 +1563,8 @@ class fields_chunk {
void initialize_with_nth_tm(int n, double kz);
// dft.cpp
void update_dfts(double timeE, double timeH, int current_step);
double dft_fields_norm2() const;
double dft_minfreq() const;

void changing_structure();
};
Expand Down Expand Up @@ -1975,6 +1979,8 @@ class fields {
dft_chunk *add_dft(const volume_list *where, const std::vector<double> &freq,
bool include_dV = true);
void update_dfts();
double dft_fields_norm();
double dft_minfreq() const;
dft_flux add_dft_flux(const volume_list *where, const double *freq, size_t Nfreq,
bool use_symmetry = true, bool centered_grid = true,
int decimation_factor = 1);
Expand Down