.. module:: magnumfe State and Materials =================== The :class:`State` class plays a central role in every simulation. It's main purpose is to keep track of the current simulation state, such as the magentization and the time. Moreover, the :code:`state` object holds the mesh along with domain information and helps to connect the different solver modules through virtual attributes. Domain Management +++++++++++++++++ In the simplest case, a :code:`state` object is generated with a plain :code:`mesh` object .. code:: python mesh = Mesh("sample.xml") state = State(mesh) Typical structures simulated with magnum.fe have a size of nanometers to microns. Since magnum.fe expects all quantities to be given in SI units, the mesh size should be given in meters. However, meshes with dimensions of :math:`1 \cdot 10^{-9}` lead to problems both at mesh generation and solution of finite-element systems. Hence it is recommended to use nanometers as mesh scale and set a proper scaling in the :class:`State` initializer. .. code:: python mesh = Mesh("sample.xml") state = State(mesh, scale = 1e-9) There are two possibilties to define domains in the mesh. Either the domain information is embedded directly in the mesh XML file or the domain markers are provided in separate :class:`MeshFunction` objects. For details, see :ref:`meshes`. Domains in the mesh files are addressed by integer IDs. The :class:`State` class can be initialized with a domain name mapping for convenient access A more convenient way to address domains is by In order to address domains by given names the :class:`State` class can be provided with a mapping of domain names to IDs. This significantly increases the readability of simulation scripts. Moreover, some domain names are of special importance since different magnum.fe module solve problems on certain regions, e.g. the LLG is by default solved on the region named 'magnetic'. Naming domains is done on intialization of the :code:`state` object. .. code:: python state = State(mesh, cell_domains = { 'magnetic': (2,4), 'conducting': (1,2,3,4,5), 'fixed' : 2, 'free' : 4 }, facet_domains = { 'contact1': 1, 'contact2': 2 }, scale = 1e-9) The naming of cell domains (volumes) and facet domains (areas) is defined by distinct Python dictionaries. The given names can be used in a variety of methods of the :class:`State` class. .. code:: python # integrate function f over 'fixed' domain f = Function(state.FunctionSpace()) assemble(f * state.dx('fixed')) # Compute volume of 'free' domain V = state.volume('free') # Compute area of 'contact1' domain A = state.area('contact1') Besides these generic uses of the domain names almost every solver class (e.g. field terms, time integrators) can be configured to solve on certain domains. For instance the :class:`LLGCvode` class solves the LLG on the 'magnetic' region be default, but can easily be configured to solve on another region. The same applies to the effective field contributions. .. code:: python # solve the LLG with exchange field only on 'some_region' exchange = ExchangeField(region = 'some_region') llg = LLGCvode([exchange], region = 'some_region') Furthermore, nonlocal field terms such as the demagnetization field and the Oersted field can be configured with a source and a target region. .. code:: python # solve the demag field generated by region 'a' on region 'b' demag = DemagField(source_region = 'a', target_region = 'b') This feature is for instance useful for the simulation of multilayer structures with a free and fixed magnetic layer that acts as spin polarizer. By restricting the source region to the free layer, the demagnetization generated by the fixed layer can be neglected while still considering its magnetic properties for the computation of the spin accumulation. When using the FEM/BEM method for the solution of the nonlocal field terms, see :ref:`open-boundary`, an additional :code:`compute_region` can be defined. This region has to contain both the source and the target region and defines the area that is solved with pure FEM. In the case of multilayer structures with alternating magnetic and nonmagnetic layers, setting the compute region to 'all' usually reduces the boundary that has to be considered for the boundary element method. Furthermore, by an appropriate choice of the compute region, small distances of boundary elements, that lead to inaccurate BEM coefficients, can be avoided. For further information refer to the reference of the :class:`DemagField` and :class:`OerstedField` class. Material Definition +++++++++++++++++++ Material parameters are managed by the :class:`State` class and made accessible via the :code:`material` attribute. They can be set either globally or in a specific domain. .. code:: python # global assignment of multiple parameters state.material = Material( alpha = 1.0, ms = 8e5 ) # global assignment of single paramter state.material.alpha = 0.1 # assign multiple parameters in region named 'fixed' state.material['fixed'] = Material( K_uni = 1e6, K_uni_axis = (1, 0, 0) ) # assign single parameter in region named 'free' state.material['free'].K_uni = 0.0 # assign spatially varying parameter state.material['magnetic'].K_uni = Expression("x[0]") Depending on the applied modules, a multitude of material parameters has to be set. Internally, all material parameters are discretized by elementwise constant functions. The parameters are applied in the order of appearance in the script, so the subsequent call of .. code:: python state.material.alpha = 1.0 state.material.alpha = 0.1 will lead to :code:`alpha` = 0.1. The attribute getters of :code:`state.material` return a FEniCS :class:`Function` object, so proper assignment can be checked by writing the material function to a file. .. code:: python File("alpha.pvd") << state.material.alpha Attributes ++++++++++ The :code:`state` object is used to communicate any simulation state between the different modules and classes of magnum.fe. Different modules require the state to provide different attributes. For instance the :class:`DemagField` module requires the :code:`state` object to provide the magnetization :code:`state.m`, while the :class:`OerstedField` module depends on the electric current :code:`state.j`. Setting an attribute on the :code:`state` object is as easy as calling :code:`state.j = Constant((1.0, 0.0, 0.0))`. Note that attributes can be either simple scalars or :class:`GenericFunction` objects, see `FEniCS manual `_. :class:`GenericFunction` is a superclass to :class:`Constant`, :class:`Expression` and :class:`Function` objects. Any :class:`GenericFunction` object will be automatically converted to a :class:`Function`. .. code:: python # Each of these commands will result in state.m being a Function object state.m = Constant((1, 0, 0)) state.m = Expression(("x[0]", "0", "0")) state.m = Function(...) :class:`GenericFunction` object are not only converted to :class:`Function` objects, but also decorated with some useful methods that are aware of the named domains. .. code:: python # Return function m restricted on the 'fixed' region state.m['fixed'] # Compute componentwise average of m state.m.average() # Compute componentwise average of m on 'free' region state.m.average('free') .. _virtual_attributes: Virtual Attributes ++++++++++++++++++ Once set, simple attributes as introduced in the preceding section remain constant. However, attributes are often required to change depending on the simulation state. For instance, in case of an AC current, the electric current obviously depends on the time. In order to provide attributes depending on other attributes, magnum.fe provides the concept of virtual attributes. Virtual attribues can be defined by a function that takes the state as single parameter and returns either a scalar or a :class:`GenericFunction`. .. code:: python # Define an AC current with amplitude A and frequency freq state.j = lambda state: Constant((A * sin(2*pi*freq*state.t), 0, 0)) The main purpose of virtual attributes is to provide certain attributes to certain solvers, e.g. the :class:`OerstedField` requires :code:`state.j` to be defined. However, you may define arbitrary virtual attributes for your own convenience. While virtual attributes are very useful to define simple time dependencies, e.g. for the definition of an AC current as shown above, they might also be used to describe complicated dependencies of the attribute from the complete simulation state. For instance the electric current might depend on a potential difference between two contacts and additionally on the magnetization configuration when dealing with spin diffusion models. magnum.fe offers a number of predefined functions that can be used as virtual attributes. Simply speaking virtual attributes provide an interface for the flexible combination of different models. .. code:: python # Compute currect j with spin diffusion model spin_acc = SpinAccumulationForPotential() state.j = spin_acc.j() # Compute Oersted field oersted = OerstedField() oersted.h(state) # Alternatively register Oersted field computation as virtual attribute state.h = oersted.h() Sometimes multiple modules will require the same virtual attribute in a single time step. If the computation of the virtual attribute is expensive, the virtual attribute should only be computed once and the result should be reused as long as possible. For this purpose virtual attribues provide a simple caching mechanism. If the virtual attributes depends on one are more other attributes, caching can be enabled by returning not the result itself, but a tuple containing the results as first entry and the names of the depending attributes as further entries. .. code:: python # compute energy -inner(m, h) and recompute only on change of state.m state.E = lambda state: (assemble(-inner(state.m, h)*dx), "m") # set m state.m = Constant((1.0, 0.0, 0.0)) energy = state.E # trigger computation of E energy = state.E # takes E from cache # reset m state.m = Constant((0.0, 1.0, 0.0)) energy = state.E # trigger computation of E Predefined virtual attributes of magnum.fe usually implement reasonable caching so you dont' have to worry about that. Masked State ++++++++++++ Complex micromagnetic models might consist of a multitude of regions and some subproblems, i.e. partial differential equations, might be defined only on certains regions. For this use case, magnum.fe provides the possibility to mask a state to a certain domain. A masked state acts exactly like a regular state, except that the mesh and all function valued attributes are region to the masked region. .. code:: python # Create a masked state for region 'magnetic' masked_state = state['magnetic'] # Override values of attribute j only in 'magnetic' region masked_state.j = Constant((1.0, 0.0, 0.0)) # This returns True assemble(Constant(1.)*state.dx('magnetic')) == assemble(Contant(1.)*masked_state.dx()) Class Reference +++++++++++++++ :class:`State` -------------- .. autoclass:: State :members: :class:`Material` ----------------- .. autoclass:: Material :members: