NanoComp / meep

free finite-difference time-domain (FDTD) software for electromagnetic simulations
GNU General Public License v2.0
1.17k stars 598 forks source link

Results From Near2Far Are Inconsistent When Loading Chunk Layout #2628

Open AwkwardTurtle42 opened 10 months ago

AwkwardTurtle42 commented 10 months ago

I've had need to run fairly lengthy meep simulations and found it more convenient to run a simulation, save the Near2Far data, and then at a later time extract whatever farfields I needed. I have also been dumping and loading the chunk layout, due to testing "split_chunks_evenly=False".

The farfield results were not what I expected, even though the near field results looked normal. I created a super simple "toy model" to make sure I wasn't missing anything obvious, a plane wave propagating through an empty 3D box with a single Near2FarRegion at the far end.

SimLayoutXZ

I still saw saw strange, asymmetric farfields and so with a bit more digging I discovered that the issues only occurred while running on multiple cores, and only while loading in the chunk_layout. It also might only occur while the simulation area is not a perfect square.

FarFieldExample

Code below:

import meep as mp
import numpy as np
import matplotlib.pyplot as plt

resolution = 4
courant = 0.5

delFreq = 0.1
centerFreq = 1/4
minFreq = centerFreq-delFreq/2
maxFreq = centerFreq+delFreq/2
max_wavelength = 1/minFreq

depth = 6
width = 6
length = 12
incident_height = max_wavelength/2
exit_height = max_wavelength/2
pml_height = max_wavelength/2

pml_layers = [mp.PML(pml_height)]

sz = pml_height + incident_height + depth + exit_height + pml_height
sx = pml_height + length + pml_height
sy = pml_height + width + pml_height

cell = mp.Volume(center=mp.Vector3(0,0,sz/2),size=mp.Vector3(sx,sy,sz))

farfield_cell = mp.Volume(center=mp.Vector3(0,0,sz-pml_height),size=(length,width,0))

geometry = []

sources = [mp.Source(mp.GaussianSource(centerFreq,fwidth=delFreq,is_integrated=True),
                component=mp.Ex,
                volume=mp.Volume(center=(0,0,pml_height),size=(sx,sy,0)))]

sim = mp.Simulation(cell_size=cell.size,
                    geometry_center=cell.center,
                    boundary_layers=pml_layers,
                    geometry=geometry,
                    sources=sources,
                    resolution=resolution,
                    Courant=courant)

n2f_obj = sim.add_near2far(centerFreq, 0, 1, mp.Near2FarRegion(volume=farfield_cell,weight=+1))

plot_slice = mp.Volume(center=(0,0,sz/2),size=(sx,0,sz))

if mp.am_master():
    plt.figure()
    sim.plot2D(output_plane=plot_slice)
    plt.savefig('SimLayoutXZ.png')

sim.run(until_after_sources=mp.stop_when_energy_decayed(dt=10, decay_by=1e-6))

sim.save_near2far("FarField", n2f_obj)
sim.dump_chunk_layout("ChunkLayout")

sim.reset_meep()

sim = mp.Simulation(cell_size=cell.size,
                    geometry_center=cell.center,
                    boundary_layers=pml_layers,
                    geometry=geometry,
                    sources=sources,
                    resolution=resolution,
                    Courant=courant,
                    # chunk_layout="ChunkLayout"
                    )

n2f_obj = sim.add_near2far(centerFreq, 0, 1, mp.Near2FarRegion(volume=farfield_cell,weight=+1))

sim.load_near2far("FarField",n2f_obj)

ff = sim.get_farfields(n2f_obj,
                        2,
                        size=mp.Vector3(sx,sy,0),
                        center=mp.Vector3(0,0,sz-pml_height+1))

np.save("Multi_Core_NoLoadChunk.npy",ff,allow_pickle=True)