tidy3d.Boundary#
- class Boundary[source]#
Bases:
Tidy3dBaseModel
Boundary conditions at the minus and plus extents along a dimension.
- Parameters:
attrs (dict = {}) – Dictionary storing arbitrary metadata for a Tidy3D object. This dictionary can be freely used by the user for storing data without affecting the operation of Tidy3D as it is not used internally. Note that, unlike regular Tidy3D fields,
attrs
are mutable. For example, the following is allowed for setting anattr
obj.attrs['foo'] = bar
. Also note that Tidy3D` will raise aTypeError
ifattrs
contain objects that can not be serialized. One can check ifattrs
are serializable by callingobj.json()
.plus (Union[Periodic, PECBoundary, PMCBoundary, PML, StablePML, Absorber, BlochBoundary] = PML(attrs={}, name=None, type='PML', num_layers=12, parameters=PMLParams(attrs={},, sigma_order=3,, sigma_min=0.0,, sigma_max=1.5,, type='PMLParams',, kappa_order=3,, kappa_min=1.0,, kappa_max=3.0,, alpha_order=1,, alpha_min=0.0,, alpha_max=0.0))) – Boundary condition on the plus side along a dimension.
minus (Union[Periodic, PECBoundary, PMCBoundary, PML, StablePML, Absorber, BlochBoundary] = PML(attrs={}, name=None, type='PML', num_layers=12, parameters=PMLParams(attrs={},, sigma_order=3,, sigma_min=0.0,, sigma_max=1.5,, type='PMLParams',, kappa_order=3,, kappa_min=1.0,, kappa_max=3.0,, alpha_order=1,, alpha_min=0.0,, alpha_max=0.0))) – Boundary condition on the minus side along a dimension.
Notes
To specify individual boundary conditions along different dimensions, instead of
BoundarySpec
, this class is used, which defines theplus
andminus
boundaries along a single dimension.Example
>>> boundary = Boundary(plus = PML(), minus = PECBoundary())
See also
BoundarySpec
Specifies boundary conditions on each side of the domain and along each dimension.
PML
A standard PML along a single dimension.
- Notebooks
Attributes
Methods
absorber
([num_layers, parameters])Adiabatic absorber boundary specification on both sides along a dimension.
bloch
(bloch_vec)Bloch boundary specification on both sides along a dimension.
bloch_from_source
(source, domain_size, axis)Bloch boundary specification on both sides along a dimension based on a given source.
bloch_on_both_sides
(values)Error if a Bloch boundary is applied on only one side.
pec
()PEC boundary specification on both sides along a dimension.
periodic
()Periodic boundary specification on both sides along a dimension.
periodic_with_pec_pmc
(values)If a PBC is specified along with PEC or PMC on the other side, manually set the PBC to PEC or PMC so that no special treatment of halos is required.
periodic_with_pml
(values)Error if PBC is specified with a PML.
pmc
()PMC boundary specification on both sides along a dimension.
pml
([num_layers, parameters])PML boundary specification on both sides along a dimension.
stable_pml
([num_layers, parameters])Stable PML boundary specification on both sides along a dimension.
Inherited Common Usage
- plus#
- minus#
- classmethod bloch_on_both_sides(values)[source]#
Error if a Bloch boundary is applied on only one side.
- classmethod periodic_with_pec_pmc(values)[source]#
If a PBC is specified along with PEC or PMC on the other side, manually set the PBC to PEC or PMC so that no special treatment of halos is required.
- classmethod periodic()[source]#
Periodic boundary specification on both sides along a dimension.
Example
>>> pbc = Boundary.periodic()
- classmethod bloch(bloch_vec)[source]#
Bloch boundary specification on both sides along a dimension.
- Parameters:
bloch_vec (complex) – Normalized component of the Bloch vector in units of 2 * pi / (size along dimension) in the background medium, along the dimension in which the boundary is specified.
Example
>>> bloch = Boundary.bloch(bloch_vec=1)
- classmethod bloch_from_source(source, domain_size, axis, medium=None)[source]#
Bloch boundary specification on both sides along a dimension based on a given source.
- Parameters:
source (Union[
GaussianBeam
,ModeSource
,PlaneWave
]) – Angled source.domain_size (float) – Size of the domain in the direction normal to the Bloch boundary
axis (int) – Axis normal to the Bloch boundary
medium (
Medium
) – Background medium associated with the Bloch vector. Default: free space.
Example
>>> from tidy3d import GaussianPulse, PlaneWave, inf >>> pulse = GaussianPulse(freq0=200e12, fwidth=20e12) >>> pw_source = PlaneWave( ... size=(inf,inf,0), source_time=pulse, direction='+', angle_theta=0.2, angle_phi=0.3) >>> bloch = Boundary.bloch_from_source(source=pw_source, domain_size=5, axis=0)
- classmethod pec()[source]#
PEC boundary specification on both sides along a dimension.
Example
>>> pec = Boundary.pec()
- classmethod pmc()[source]#
PMC boundary specification on both sides along a dimension.
Example
>>> pmc = Boundary.pmc()
- classmethod pml(num_layers=12, parameters=PMLParams(attrs={}, sigma_order=3, sigma_min=0.0, sigma_max=1.5, type='PMLParams', kappa_order=3, kappa_min=1.0, kappa_max=3.0, alpha_order=1, alpha_min=0.0, alpha_max=0.0))[source]#
PML boundary specification on both sides along a dimension.
- Parameters:
num_layers (int = 12) – Number of layers of standard PML to add to + and - boundaries.
parameters (
PMLParams
) – Parameters of the complex frequency-shifted absorption poles.
Example
>>> pml = Boundary.pml(num_layers=20)
- classmethod stable_pml(num_layers=40, parameters=PMLParams(attrs={}, sigma_order=3, sigma_min=0.0, sigma_max=1.0, type='PMLParams', kappa_order=3, kappa_min=1.0, kappa_max=5.0, alpha_order=1, alpha_min=0.0, alpha_max=0.9))[source]#
Stable PML boundary specification on both sides along a dimension.
- Parameters:
num_layers (int = 40) – Number of layers of ‘stable’ PML to add to + and - boundaries.
parameters (
PMLParams
) – ‘Stable’ parameters of the complex frequency-shifted absorption poles.
Example
>>> stable_pml = Boundary.stable_pml(num_layers=40)
- classmethod absorber(num_layers=40, parameters=AbsorberParams(attrs={}, sigma_order=3, sigma_min=0.0, sigma_max=6.4, type='AbsorberParams'))[source]#
Adiabatic absorber boundary specification on both sides along a dimension.
- Parameters:
num_layers (int = 40) – Number of layers of absorber to add to + and - boundaries.
parameters (
PMLParams
) – Adiabatic absorber parameters.
Example
>>> absorber = Boundary.absorber(num_layers=40)
- __hash__()#
Hash method.