emopt.fdfd¶
This module defines a standard interface for defining finite difference frequency domain (FDFD) solvers. The FDFD method formulates the frequency domain Maxwell’s equations explicitly in the form \(A x = b\). The matrix \(A\) is built up by discretizing the materials of the system on a rectangular grid. The curl operators in Maxwell’s equations are then approximated on this grid using centered finite differences. In order to accomodate these centered differences, a set of dislocated grids (i.e. an array of Yee cells) is used. The fields are solved for on this same grid.
In order to run an FDFD simulation, the materials and sources of the system must be defined. Once these are defined, the system matrix (\(A\)) can be constructed and a simulation performed. In addition to methods for obtaining the resulting fields and calculating the source power, the FDFD interface also defines a function which exposes the diagonal elements of \(A\). This is useful for performing sensitivity analysis/calculating gradients which is a core goal of EMOpt.
In addition to the FDFD interface, a variety of specific implementations are defined. In particular:
FDFD_TE : A 2D FDFD solver with TE-polarized fields. The system is infinitely extruded in z and waves can propagate alogn x and y. The allowed sources of the system are \(J_z\), \(M_x\), \(M_y\) and the non-zeros field components are \(E_z\), \(H_x\), and \(H_y\).
FDFD_TM : A 2D FDFD solver with TM-polarized fields. The system is infinitely extruded in z and waves can propagate alogn x and y. The allowed sources of the system are \(M_z\), \(J_x\), \(J_y\) and the non-zeros field components are \(H_z\), \(E_x\), and \(E_y\).
FDFD_3D : A full vector 3D FDFD solver. All current source and field components are used. Material distributions in this solver are currently restricted to planar structures.
It is important to note that the terminology used here for ‘TE’ and ‘TM’ is consistent with the more classical electromagnetic definition and not with the silicon photonics definition.
Examples
See emopt/examples/ for detailed examples.
Todo
1. Reimplement get_field and get_field_interp in parallel and add support for Domains in the form FDFD.get_field(self, component, domain=None). The parallel support will be needed when 3D solvers are implemented regardless, so it is a good idea to implement the functionality in the already working 2D solvers. Note: this will require modification of the AdjointMethod class as well; to make this as simple as possible, the get_field function should return a correctly-sized zeroed numpy array on the non-master nodes.
2. Reimplement the set_source function with Domain support. This will be necessary in 3D when the total grid size is very large.
-
class
emopt.fdfd.FDFD(ndims)¶ Bases:
emopt.simulation.MaxwellSolverFinite difference frequency domain solver.
This class provides the generalized interface for a finite diffrence frequency domain forward and adjoint solver. It may be extended to implement 2D semi-vector and 3D full-vector Maxwell solvers.
Notes
Finite difference frequency domain solvers must express the discretized Maxwell’s equations in the form \(A x = b\). In EMOpt, solving this system of equations is handled using the petsc4py library, although this is not strictly necessary.
The majority of methods in this class are abstract and need to be implemented by a child class. The exception to this is the function
get_A_diag()which should not be necessary to override.In general, there are two different ways to access the results of an FDFD simulation. The first way is through direct calls to the
get_field()andget_field_interp()functions. The second is by specifying field domains before running a forward simulation. If field domains are provided, the fields in these domains will be saved immediately after running the simulation. This is useful in situations where the same fields will be retrieved multiple times between simulations and in cases where the fields need to be accessed from the master node without involving other higher-rank nodes.-
solve_forward(self)¶ Simulate Maxwell’s equations
-
solve_adjoint(self)¶ Simulate the transposed Maxwell’s equations
-
get_field(self, component)¶ Get the desired component of the raw (uninterpolated) field
-
get_field_interp(self, component)¶ Get the desired component of the interpolated field
-
get_adjoint_field(self, component)¶ Get the desired component of the raw (uninterpolated) adjoint field
-
get_adjoint_field_interp(self, component)¶ Get the desired component of the interpolated adjoint field
-
build(self)¶ Build the system matrix \(A\).
-
update(self)¶ Update the system matrix \(A\)
-
set_sources(self, src)¶ Set the discretized current sources used in the forward simulation.
-
set_adjoint_sources(self, src)¶ Set the discretized current sources used in the adjoint simulation.
-
get_source_power(self)¶ Get the total power generated by the sources.
-
get_A_diag(self, vdiag=None)¶ Retrieve the diagonal elements of \(A\)
-
nunks¶ The number of unknowns solved for in the simulation.
-
field_domains¶ The list of DomainCoordinates in which fields are recorded immediately following a forward solve.
-
saved_fields¶ The list of fields saved in in the stored DomainCoordinates
-
source_power¶ The source power injected into the system.
-
calc_ydAx(Adiag0)¶ Calculate y^T * (A1-A0) * x.
Parameters: Adiag0 (PETSc.Vec) – The diagonal of the FDFD matrix. Returns: The product y^T * (A1-A0) * x Return type: complex
-
get_A_diag(vdiag=None) Get the diagonal entries of the system matrix A.
Todo
Implement this using a sub-class-dependent in-place copy to speed things up?
Parameters: vdiag (petsc4py.PETSc.Vec) – Vector with dimensions Mx1 where M is equal to the number of diagonal entries in A. Returns: Return type: **(Master node only)** the diagonal entries of A.
-
nunks
-
spy_A(**kwargs)¶
-
-
class
emopt.fdfd.FDFD_3D(X, Y, Z, dx, dy, dz, wavelength, mglevels=None, rtol=1e-06, low_memory=False, verbose=True)¶ Bases:
emopt.fdfd.FDFDSimulate Maxwell’s equations in 3D.
Notes
1. Units used for length parameters do not matter as long as they are consistent with one another
2. The dimensions of the system will be modified such that they are consistent with the provided grid spacing. You should always reinitialize your scripts X, Y, Z variables with sim.X, sim.Y, and sim.Z in order to make sure everything is consistent.
3. The 3D solver internally uses an iterative Krylov subspace method with a multigrid preconditioner. This is quite a bit more complicated than the direct solver used in 2D so careful attention should be paid when setting up the 3D solver. Additional documentation on this will come in the future.
4. The solver supports both real and complex material values, however the source power calculation currently only accounts for real-valued materials. This will be updated in the future.
Parameters: - W (float) – Approximate width of the simulation region. If this is not an integer multiple of dx, this will be increased slightly
- H (float) – Approximate height of the simulation region. If this is not an integer multiple of dy, this will be increased slightly
- dx (float) – Grid spacing in the x direction.
- dy (float) – Grid spacing in the y direction.
- wavelength (float) – Vacuum wavelength of EM fields.
- mglevels (int) – The number of levels to use in the multigrid prevonditioner. If running very coarse resolutions, this value should be decreased (to 2). If running a very high resolution simulation, this value should be increased (to 4 or more depending on resolution).
- rtol (float) – The relative tolerance used to terminate the internal iterative solver. Smaller number = higher accuracy solution.
- low_memory (boolean) – If True, reduce the memory footprint of the solver at the cost of performance. This is achieved by using an iterative method for the coarse solver instead of a direct LU factorization. (default = False)
-
pml_power¶ The power of the polynomial used to define the spatial-dependence of pml_sigma.
Type: float
-
field_domains¶ The list of DomainCoordinates in which fields are recorded immediately following a forward solve.
Type: list of DomainCoordinates
-
saved_fields¶ The list of (Ez, Hx, Hy) fields saved in in the stored DomainCoordinates
Type: list of numpy.ndarray
-
w_pml¶ List of PML widths in real spatial units. This variable can be reinitialized in order to change the PML widths. The list is ordered as follows: [xmin, xmax, ymin, ymax, zmin, zmax]
Type: list of float
-
real_mats¶ If True, assume all permittivity and permeability values are real valued (in power calculations).
Type: boolean
-
Nx¶
-
Ny¶
-
Nz¶
-
X¶
-
Y¶
-
Z¶
-
bc¶
-
build()¶ (Re)Build the system matrix.
Maxwell’s equations are solved by compiling them into a linear system of the form \(Ax=b\). Here, we build up the structure of A which contains the curl operators, mateiral parameters, and boundary conditions.
This function must be called at least once after the material distributions of the system have been defined (through a call to
set_materials()).Raises: Exception– If the material distributions of the simulation have not been set prior to calling this function.
-
buildA(l)¶
-
buildRst(l)¶
-
calc_ydAx(Adiag0)¶ Calculate y^T * (A1-A0) * x.
Parameters: Adiag0 (PETSc.Vec) – The diagonal of the FDFD matrix. Returns: The product y^T * (A1-A0) * x Return type: complex
-
dx
-
dy
-
dz¶
-
get_adjoint_field(component, domain=None, squeeze=False)¶ Get the adjoint field.
Notes
1. This function is primarily intended as a diagnostics tool. If implementing and adjoint method, consult the emopt.adjoint_method documentation.
2. Currently aggregation from the different nodes is done is a slow way. It is best to avoid making many repetative calls to this function.
Parameters: - component (str) – The desired adjoint field component to retrieve (Ex, Ey, Ez, Hx, Hy, Hz). Note that the adjoint field is somewhat fictitious, so “field component” refers to the part of the solution vector from which to collect the fields.
- domain (emopt.misc.DomainCoordinates (optional)) – On field data from the specified domain is retrieved. If no domain is specified, the whole simulation domain will be used. (default=None)
Returns: The desired fields contained within the specified domain.
Return type: numpy.ndarray
-
get_field(component, domain=None, squeeze=False)¶ Get the uninterpolated field.
Notes
Currently aggregation from the different nodes is done is a slow way. It is best to avoid making many repetative calls to this function.
Parameters: - component (str) – The desired field component to retrieve (Ex, Ey, Ez, Hx, Hy, Hz)
- domain (emopt.misc.DomainCoordinates (optional)) – On field data from the specified domain is retrieved. If no domain is specified, the whole simulation domain will be used. (default=None)
Returns: The desired fields contained within the specified domain.
Return type: numpy.ndarray
-
get_field_interp(component, domain=None, squeeze=False)¶ Get the desired field component.
Internally, fields are solved on a staggered grid. In most cases, it is desirable to know all of the field components at the same sets of positions. This requires that we interpolate the fields onto a single grid. In emopt, we interpolate all field components onto the Ez grid.
Parameters: - component (str) – The desired field component.
- domain (misc.DomainCoordinates (optional)) – The domain from which the field is retrieved. (default = None)
Returns: The interpolated field
Return type: numpy.ndarray
-
get_source_power()¶ Get source power.
The source power is the total electromagnetic power radiated by the electric and magnetic current sources.
Returns: The source power. Return type: float
-
set_adjoint_sources(dFdxs)¶ Set the adjoint sources.
This function exists to maintain compatibility with the adjoint method class which expects a function which accepts a single argument in order to set the adjoint source distribution.
This functionality is complicated in 3D as we have to be a bit more careful with how we share data between processors.
Notes
- This function clears the current adjoint source vector when called.
Parameters: dFdxs (tuple of lists) – Tuple containing a list of numpy.ndarray as the first element and a list of correpsonding DomainCoordinates as the second argument.
-
set_materials(eps, mu)¶ Set the material distributions of the system to be simulated.
Parameters: - eps (emopt.grid.Material) – The spatially-dependent permittivity of the system
- mu (emopt.grid.Material) – The spatially-dependent permeability of the system
-
set_sources(src, domain, mode_num=0)¶ Set the sources of the system used in the forward solve.
Notes
1. Like the underlying fields, the current sources are represented on a set of shifted grids. If manually setting the source, this should be taken into account.
2. If using a mode source to set the sources, the mode source must be built and run prior to calling this function.
3. This function can be called repeatedly to modify the source distribution. Calling this function does not zero the current distribution.
Parameters: src (modes.ModeFullVector or tuple of numpy.ndarray) – Either a ModeFullVector object or a set of arrays containing current source distributions.
-
solve_adjoint()¶ Solve the adjoint simulation.
This is equivalent to a solution of the transpose of Maxwell’s equations, which is achieved by either iteratively or directly solving \(Ax=b\).
-
solve_forward()¶ Solve the forward simulation.
This is equivalent to a solution of Maxwell’s equations, which is achieved by either iteratively or directly solving \(Ax=b\).
After the completing the forward solve, the full solution contained in x is gathered on the master node so that the fields can be accessed.
-
update(bbox=None)¶ (Partially) Update the system matrix.
Updating the system matrix involves recomputing the material distribution of the system.
If a bounding box is provided, only that portion of the grid will be updates.
Notes
This only updates the high resolution grid.
-
update_adjoint_sources(src, domain)¶ Update the adjoint source.
Parameters: - src (numpy.ndarray) – The source distribution
- domain (misc.DomainCoordinates) – The region in the grid corresponding to src.
-
update_multigrid(l)¶ Update the material distributions of the restricted system matrices
Notes
- l : int
- The multigrid level to update.
-
update_saved_fields()¶ Update the saved fields based on the domains contained in FDFD_3D.field_domains.
Notes
- This is primarily for internal use only.
- The fields are NOT squeezed
-
w_pml
-
class
emopt.fdfd.FDFD_TE(X, Y, dx, dy, wavelength, solver='auto', ksp_solver='gmres')¶ Bases:
emopt.fdfd.FDFDSimulate Maxwell’s equations in 2D with TE-polarized fields.
Notes
1. Units used for length parameters do not matter as long as they are consistent with one another
2. The width and the height of the system will be modified in order to ensure that they are an integer multiple of dx and dy. It is important that this modified width and height be used in any future calculations.
Parameters: - X (float) – Approximate width of the simulation region. If this is not an integer multiple of dx, this will be increased slightly
- Y (float) – Approximate height of the simulation region. If this is not an integer multiple of dy, this will be increased slightly
- dx (float) – Grid spacing in the x direction.
- dy (float) – Grid spacing in the y direction.
- wavelength (float) – Vacuum wavelength of EM fields.
- solver (str) – The type of solver to use. The possible options are ‘direct’ (direct LU solver), ‘iterative’ (unpreconditioned iterative solver), ‘iterative_lu’ (iterative solver with LU preconditioner, or ‘auto’. (default = ‘auto’)
- ksp_solver (str) – The type of Krylov subspace solver to use. See the petsc4py documentation for the possible types. Note: this flag is ignored if ‘direct’ or ‘auto’ is chosen for the solver parameter. (default = ‘gmres’)
- verbose (boolean (Optional)) – Indicates whether or not progress info should be printed on the master node. (default=True)
-
field_domains¶ The list of DomainCoordinates in which fields are recorded immediately following a forward solve.
Type: list of DomainCoordinates
-
saved_fields¶ The list of (Ez, Hx, Hy) fields saved in in the stored DomainCoordinates
Type: list of numpy.ndarray
-
w_pml¶ List of PML widths in real spatial units. This variable can be reinitialized in order to change the PML widths
Type: list of float
-
real_materials¶ If True, assume that epsilon and mu are real-valued. This can speed up some calculations. (defaults to False)
Type: bool
-
H The height of the systme. This height will only exactly match the height passed during initialization if it is equal to an integer multiple of grid cells.
-
Hreal¶ Height of the simulation, excluding PMLs.
-
M
-
N
-
W The width of the systme. This width will only exactly match the width passed during initialization if it is equal to an integer multiple of grid cells.
-
Wreal¶ Width of the simulation, excluding PMLs.
-
X This is needed for MaxwellSolver implementation.
Todo
Replace self._W with self._X
-
Xreal¶ Width of the simulation, excluding PMLs.
-
Y Also needed for MaxwellSolver implementation.
Todo
Replace self._H with self._Y
-
Yreal¶ Height of the simulation, excluding PMLs.
-
bc¶
-
build()¶ (Re)Build the system matrix.
Maxwell’s equations are solved by compiling them into a linear system of the form \(Ax=b\). Here, we build up the structure of A which contains the curl operators, mateiral parameters, and boundary conditions.
This function must be called at least once after the material distributions of the system have been defined (through a call to
set_materials()).Notes
1. In the current Yee cell configuration, the \(E_z\)’s are positioned at the center of the cell, the \(H_x\)’s are shifted in the positive y direction by half of a grid cell, and the \(H_y\)’s are shifted in the negative x direction by half a grid cell.
2. This function can be rather slow for larger problems. In general, callin this function more than once should be avoided. Instead, the
update()function should be sufficient when only the material distributions of the system have been modified3. The current PML implementation is not likely ideal. Nonetheless, it seems to work ok.
Raises: Exception– If the material distributions of the simulation have not been set prior to calling this function.
-
calc_ydAx(Adiag0)¶ Calculate y^T * (A1-A0) * x.
Parameters: Adiag0 (PETSc.Vec) – The diagonal of the FDFD matrix. Returns: The product y^T * (A1-A0) * x Return type: complex
-
dx
-
dy
-
eps¶
-
get_adjoint_field(component, domain=None)¶ Get the desired raw adjoint field component.
Notes
this function only returns a non-empty field on the master node.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Ez’, ‘Hx’, or ‘Hy’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Ez’, ‘Hx’, or ‘Hy’Returns: (Master node only) A numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_field(component, domain=None)¶ Get the desired (forward) field component.
This function returns the RAW field which is defined on the dislocated grids. In most situations, the
get_field_interp()function should be prefered.Notes
This function only returns a non-empty field on the master node.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Ez’, ‘Hx’, or ‘Hy’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Ez’, ‘Hx’, or ‘Hy’Returns: (Master node only) A 2D numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_field_interp(component, domain=None)¶ Get the desired (forward) interpolated field.
When solving Maxwell’s equations on a rectangular grid, we actually set up a set of dislocated grids in which the electric and magnetic field components are computed at slightly different positions in space. As as result, it is often convenient or even necessary to interpolate the fields such that they are all known at the same points in space. This is particularly important when computing power-related quantities.
For the sake of simplicity, a simple linear interpolation scheme is used to compute Hx and Hy at the position of Ez. Ez is the same as the uninterpolated Ez field.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Ez’, ‘Hx’, or ‘Hy’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Ez’, ‘Hx’, or ‘Hy’Returns: (Master node only) A 2D numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_source_power()¶ Get the source power.
In general, FDFD_TE.source_power should be used to access the source power instead. The source power is computed following every forward solve, so calling this function explicitly is generally unnecessary.
Notes
- The source power is computed using the interpolated fields.
Returns: Electromagnetic power generated by source. Return type: float
-
mu¶
-
set_adjoint_sources(src)¶ Set the sources of the system used in the adjoint solve.
The details of these sources are identical to the forward solution current sources.
Parameters: src (tuple of numpy.ndarray) – The current sources in the form (Jz_adj, Mx_adj, My_adj). Each array in the tiple should be a 2D numpy.ndarry with dimensions MxN.
-
set_materials(eps, mu)¶ Set the material distributions of the system to be simulated.
Parameters: - eps (emopt.grid.Material) – The spatially-dependent permittivity of the system
- mu (emopt.grid.Material) – The spatially-dependent permeability of the system
-
set_sources(src, src_domain=None, mindex=0)¶ Set the sources of the system used in the forward solve.
The sources can be defined either using three numpy arrays or using a mode source. When using a mode source, the corresponding current density arrays will be automatically generated and extracted.
Notes
Like the underlying fields, the current sources are represented on a set of shifted grids. In particular, \(J_z\)’s are all located at the center of a grid cell, the \(M_x\)’s are shifted in the positive y direction by half a grid cell, and the \(M_y\)’s are shifted in the negative x direction by half of a grid cell. It is important to take this into account when defining the current sources.
Todo
1. Implement a more user-friendly version of these sources (so that you do not need to deal with the Yee cell implementation).
- Implement this in a better parallelized way
Parameters: - src (tuple of numpy.ndarray or modes.ModeTE) – (Option 1) The current sources in the form (Jz, Mx, My). Each array in the tuple should be a 2D numpy.ndarry with dimensions MxN. (Option 2), a mode source which has been built and solved.
- src_domain (emopt.misc.DomainCoordinates (optional)) – Specifies the location of the provided current source distribution or mode source. If None, it is assumed that source arrays have been provided and those source arrays span the whole simulation region (i.e. have size MxN)
- mindx (int (optional)) – Specifies the index of the mode source (only used if a mode source is passed in as src)
-
solve_adjoint()¶ Solve the adjoint simulation.
The adjoint simulation a transposed version of the forward simulation of the form \(A^T u = c\). The adjoint simulation is very useful for calculating sensitivies of the structure using the adjoint method.
-
solve_forward()¶ Solve the forward simulation.
This is equivalent to a solution of Maxwell’s equations, which is achieved by either iteratively or directly solving \(Ax=b\).
After the completing the forward solve, the full solution contained in x is gathered on the master node so that the fields can be accessed.
-
solver_type¶
-
update(bbox=None)¶ Update only the material values stored in A.
In many situations, it is desirable to be able to modify the physical structure of the simulation. Doing so only requires updating the material data in A, which corresponds to the diagonal elements. This is a much faster.
In many cases, only a small part of the full underlying grid needs to be updated. In order to handle this, an optional bounding box can be passed to limit the extent of the update.
To further improve performance, the bulk of this function has been implemented directly in C++. Refer to Grid_ctypes.hpp for details.
Parameters: bbox (tuple of ints) – bounding box defining the block of x and y indices to update in the format (x1, x2, y1, y2). If None, then the whole grid is updated. (default = None)
-
update_saved_fields()¶
-
w_pml
-
w_pml_bottom
-
w_pml_left
-
w_pml_right
-
w_pml_top
-
wavelength
-
class
emopt.fdfd.FDFD_TM(X, Y, dx, dy, wavelength, solver='auto', ksp_solver='gmres')¶ Bases:
emopt.fdfd.FDFD_TESimulate Maxwell’s equations in 2D with TM-polarized fields.
Notes
1. Units used for length parameters do not matter as long as they are consistent with one another
2. The width and the height of the system will be modified in order to ensure that they are an integer multiple of dx and dy. It is important that this modified width and height be used in any future calculations.
Parameters: - X (float) – Approximate width of the simulation region. If this is not an integer multiple of dx, this will be increased slightly
- Y (float) – Approximate height of the simulation region. If this is not an integer multiple of dy, this will be increased slightly
- dx (float) – Grid spacing in the x direction.
- dy (float) – Grid spacing in the y direction.
- wavelength (float) – Vacuum wavelength of EM fields.
- solver (str) – The type of solver to use. The possible options are ‘direct’ (direct LU solver), ‘iterative’ (unpreconditioned iterative solver), ‘iterative_lu’ (iterative solver with LU preconditioner, or ‘auto’. (default = ‘auto’)
- ksp_solver (str) – The type of Krylov subspace solver to use. See the petsc4py documentation for the possible types. Note: this flag is ignored if ‘direct’ or ‘auto’ is chosen for the solver parameter. (default = ‘gmres’)
- verbose (boolean (Optional)) – Indicates whether or not progress info should be printed on the master node. (default=True)
-
field_domains¶ The list of DomainCoordinates in which fields are recorded immediately following a forward solve.
Type: list of DomainCoordinates
-
saved_fields¶ The list of (Ez, Hx, Hy) fields saved in in the stored DomainCoordinates
Type: list of numpy.ndarray
-
w_pml¶ List of PML widths in real spatial units. This variable can be reinitialized in order to change the PML widths
Type: list of float
-
bc¶
-
build()¶ (Re)Build the system matrix.
Maxwell’s equations are solved by compiling them into a linear system of the form \(Ax=b\). Here, we build up the structure of A which contains the curl operators, mateiral parameters, and boundary conditions.
This function must be called at least once after the material distributions of the system have been defined (through a call to
set_materials()).Notes
1. In the current Yee cell configuration, the \(H_z\)’s are positioned at the center of the cell, the \(E_x\)’s are shifted in the positive y direction by half of a grid cell, and the \(E_y\)’s are shifted in the negative x direction by half a grid cell.
2. This function can be rather slow for larger problems. In general, callin this function more than once should be avoided. Instead, the
update()function should be sufficient when only the material distributions of the system have been modified3. The current PML implementation is not likely ideal. Nonetheless, it seems to work ok.
Raises: Exception– If the material distributions of the simulation have not been set prior to calling this function.
-
eps¶
-
get_A_diag(vdiag=None)¶ Get the diagonal entries of the system matrix A.
Parameters: vdiag (petsc4py.PETSc.Vec) – Vector with dimensions Mx1 where M is equal to the number of diagonal entries in A. Returns: Return type: **(Master node only)** the diagonal entries of A.
-
get_adjoint_field(component, domain=None)¶ Get the desired raw adjoint field component.
Notes
this function only returns a non-empty field on the master node.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Hz’, ‘Ex’, or ‘Ey’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Hz’, ‘Ex’, or ‘Ey’Returns: (Master node only) A numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_field(component, domain=None)¶ Get the desired (forward) field component.
This function returns the RAW field which is defined on the dislocated grids. In most situations, the
get_field_interp()function should be prefered.Notes
This function only returns a non-empty field on the master node.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Hz’, ‘Ex’, or ‘Ey’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Hz’, ‘Ex’, or ‘Ey’Returns: (Master node only) A 2D numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_field_interp(component, domain=None)¶ Get the desired (forward) interpolated field.
When solving Maxwell’s equations on a rectangular grid, we actually set up a set of dislocated grids in which the electric and magnetic field components are computed at slightly different positions in space. As as result, it is often convenient or even necessary to interpolate the fields such that they are all known at the same points in space. This is particularly important when computing power-related quantities.
For the sake of simplicity, a simple linear interpolation scheme is used to compute Ex and Ey at the position of Hz. Hz is the same as the uninterpolated Hz field.
Parameters: - component (str) – The desired field component to be retrieved (either ‘Hz’, ‘Ex’, or ‘Ey’)
- domain (emopt.misc.DomainCoordinates) – The domain in which the field is retrieved
Raises: ValueError– If the supplied component is not ‘Hz’, ‘Ex’, or ‘Ey’Returns: (Master node only) A 2D numpy.ndarray containing the desired field component
Return type: numpy.ndarray
-
get_source_power()¶ Get the source power.
Notes
- The source power is computed using the interpolated fields.
Returns: Electromagnetic power generated by source. Return type: float
-
mu¶
-
set_adjoint_sources(src)¶ Set the sources of the system used in the adjoint solve.
The details of these sources are identical to the forward solution current sources.
Parameters: src (tuple of numpy.ndarray) – The current sources in the form (Mz_adj, Jx_adj, Jy_adj). Each array in the tiple should be a 2D numpy.ndarry with dimensions MxN.
-
set_materials(eps, mu)¶ Set the material distributions of the system to be simulated.
Parameters: - eps (emopt.grid.Material) – The spatially-dependent permittivity of the system
- mu (emopt.grid.Material) – The spatially-dependent permeability of the system
-
set_sources(src, src_domain=None, mindex=0)¶ Set the sources of the system used in the forward solve.
The sources can be defined either using three numpy arrays or using a mode source. When using a mode source, the corresponding current density arrays will be automatically generated and extracted.
Notes
Like the underlying fields, the current sources are represented on a set of shifted grids. In particular, \(M_z\)’s are all located at the center of a grid cell, the \(J_x\)’s are shifted in the positive y direction by half a grid cell, and the \(J_y\)’s are shifted in the negative x direction by half of a grid cell. It is important to take this into account when defining the current sources.
Todo
1. Implement a more user-friendly version of these sources (so that you do not need to deal with the Yee cell implementation).
- Implement this in a better parallelized way
Parameters: - src (tuple of numpy.ndarray or modes.ModeTM) – (Option 1) The current sources in the form (Mz, Jx, Jy). Each array in the tuple should be a 2D numpy.ndarry with dimensions MxN. (Option 2), a mode source which has been built and solved.
- src_domain (emopt.misc.DomainCoordinates (optional)) – Specifies the location of the provided current source distribution or mode source. If None, it is assumed that source arrays have been provided and those source arrays span the whole simulation region (i.e. have size MxN)
- mindx (int (optional)) – Specifies the index of the mode source (only used if a mode source is passed in as src)
-
update_saved_fields()¶