|
TDMS
Time Domain Maxwell Solver
|
#include <field.h>
Public Member Functions | |
| Field ()=default | |
| Field (int I_total, int J_total, int K_total) | |
| void | allocate () |
| void | zero () |
| void | allocate_and_zero () |
| void | normalise_volume () |
| Normalises the field entries by dividing by the angular norm. | |
| void | set_phasors (SplitField &F, int n, double omega, double dt, int Nt) |
| void | add_to_angular_norm (int n, int Nt, SimulationParameters ¶ms) |
| Compute the phasor_norm of the current field and add it to the current norm-value. | |
| std::complex< double > | phasor_norm (double f, int n, double omega, double dt, int Nt) |
| virtual double | phase (int n, double omega, double dt)=0 |
| virtual std::complex< double > | interpolate_to_centre_of (AxialDirection d, CellCoordinate cell)=0 |
| Interpolates a Field component to the centre of a Yee cell. | |
| void | interpolate_over_range (mxArray *x_out, mxArray *y_out, mxArray *z_out, int i_lower, int i_upper, int j_lower, int j_upper, int k_lower, int k_upper, Dimension mode=Dimension::THREE) |
| Interpolates the Field over the range provided. | |
| void | interpolate_over_range (mxArray *x_out, mxArray *y_out, mxArray *z_out, Dimension mode=Dimension::THREE) |
| Interpolates the Field over the range provided. | |
| virtual void | interpolate_transverse_electric_components (CellCoordinate cell, std::complex< double > *x_at_centre, std::complex< double > *y_at_centre, std::complex< double > *z_at_centre)=0 |
| Interpolates the Field's transverse electric components to the centre of Yee cell i,j,k. | |
| virtual void | interpolate_transverse_magnetic_components (CellCoordinate cell, std::complex< double > *x_at_centre, std::complex< double > *y_at_centre, std::complex< double > *z_at_centre)=0 |
| Interpolates the Field's transverse magnetic components to the centre of Yee cell i,j,k. | |
| void | set_values_from (Field &other) |
| double | normalised_difference (Field &other) |
| Computes the maximum pointwise absolute difference of the other field to this one, divided by the largest absolute value of this field's components. | |
| Public Member Functions inherited from Grid | |
| int | max_IJK_tot () const |
| void | set_preferred_interpolation_methods (tdms_flags::InterpolationMethod im) |
| Set the preferred interpolation methods. | |
Public Attributes | |
| double | ft = 0. |
| std::complex< double > | angular_norm = 0. |
| XYZTensor3D< double > | real |
| XYZTensor3D< double > | imag |
| int | il = 0 |
| int | iu = 0 |
| int | jl = 0 |
| int | ju = 0 |
| int | kl = 0 |
| int | ku = 0 |
| Public Attributes inherited from Grid | |
| IJKDimensions | tot = {0, 0, 0} |
Additional Inherited Members | |
| Protected Attributes inherited from Grid | |
| tdms_flags::InterpolationMethod | interpolation_method |
A complex field defined over a grid. Has real and imaginary components at each (x, y, z) grid point
|
default |
Default no arguments constructor
| Field::Field | ( | int | I_total, |
| int | J_total, | ||
| int | K_total ) |
Constructor of the field with a defined size in the x, y, z Cartesian dimensions
| void Field::add_to_angular_norm | ( | int | n, |
| int | Nt, | ||
| SimulationParameters & | params ) |
Compute the phasor_norm of the current field and add it to the current norm-value.
Note that the calls: add_to_angular_norm(n, Nt, params), and add_to_angular_norm(n % Nt, Nt, params),
are equivalent up to numerical precision when using n = tind, and Nt = Nsteps as set in the main simulation.
| n | The number of timesteps INTO the current sinusoidal period |
| Nt | The number of timesteps in a sinusoidal period |
| params | The simulation parameters |
| void Field::allocate | ( | ) |
Allocate the memory appropriate for all the 3D tensors associated with this split field
|
inline |
Allocate and set to zero all components of the field
| void Field::interpolate_over_range | ( | mxArray * | x_out, |
| mxArray * | y_out, | ||
| mxArray * | z_out, | ||
| Dimension | mode = Dimension::THREE ) |
Interpolates the Field over the range provided.
Default range is to interpolate to the midpoint of all consecutive points.
| [out] | x_out,y_out,z_out | Output arrays for interpolated values |
| i_lower,j_lower,k_lower | Lower index for interpolation in the i,j,k directions, respectively | |
| i_upper,j_upper,k_upper | Upper index for interpolation in the i,j,k directions, respectively | |
| mode | Determines which field components to compute, based on the simulation Dimension |
| void Field::interpolate_over_range | ( | mxArray * | x_out, |
| mxArray * | y_out, | ||
| mxArray * | z_out, | ||
| int | i_lower, | ||
| int | i_upper, | ||
| int | j_lower, | ||
| int | j_upper, | ||
| int | k_lower, | ||
| int | k_upper, | ||
| Dimension | mode = Dimension::THREE ) |
Interpolates the Field over the range provided.
Default range is to interpolate to the midpoint of all consecutive points.
| [out] | x_out,y_out,z_out | Output arrays for interpolated values |
| i_lower,j_lower,k_lower | Lower index for interpolation in the i,j,k directions, respectively | |
| i_upper,j_upper,k_upper | Upper index for interpolation in the i,j,k directions, respectively | |
| mode | Determines which field components to compute, based on the simulation Dimension |
|
pure virtual |
Interpolates a Field component to the centre of a Yee cell.
| d | Field component to interpolate |
| cell | Index (i,j,k) of the Yee cell to interpolate to the centre of |
Implemented in ElectricField, and MagneticField.
|
pure virtual |
Interpolates the Field's transverse electric components to the centre of Yee cell i,j,k.
| [in] | cell | Yee cell index |
| [out] | x_at_centre,y_at_centre,z_at_centre | Addresses to write interpolated values for the x,y,z components (respectively) |
Implemented in ElectricField, and MagneticField.
|
pure virtual |
Interpolates the Field's transverse magnetic components to the centre of Yee cell i,j,k.
| [in] | cell | Yee cell index |
| [out] | x_at_centre,y_at_centre,z_at_centre | Addresses to write interpolated values for the x,y,z components (respectively) |
Implemented in ElectricField, and MagneticField.
| void Field::normalise_volume | ( | ) |
Normalises the field entries by dividing by the angular norm.
Specifically, real[c][k][j][i] + i imag[c][k][j][i] = ( real[c][k][j][i] + i imag[c][k][j][i] ) / angular_norm
| double Field::normalised_difference | ( | Field & | other | ) |
Computes the maximum pointwise absolute difference of the other field to this one, divided by the largest absolute value of this field's components.
Specifically, we compute and return the quantity ![$ \frac{ \max_{i,j,k}(this[k][j][i] - other[k][j][i]) }{
\max_{i,j,k}this[k][j][i] } $](form_0.png)
The other field must have the same dimensions as this field.
| other | The other field |
| void Field::set_phasors | ( | SplitField & | F, |
| int | n, | ||
| double | omega, | ||
| double | dt, | ||
| int | Nt ) |
Set the phasors for this field, given a split field. Result gives field according to the exp(-iwt) convention
| F | The split-field to read values from |
| n | The current timestep |
| omega | Angular frequency |
| dt | Timestep |
| Nt | Number of timesteps in a sinusoidal period |
| void Field::set_values_from | ( | Field & | other | ) |
Set the values of all components in this field from another, equally sized field
| void Field::zero | ( | ) |
Set all the values of all components of the field to zero
| int Field::il = 0 |
Upper (u) and lower (l) indices in the x,y,z directions. e.g. il is the first non-pml cell in the i direction and iu the last in the corresponding split field grid