forked from lijiext/lammps
120 lines
8.7 KiB
HTML
120 lines
8.7 KiB
HTML
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN">
|
|
<html><head><meta http-equiv="Content-Type" content="text/html;charset=iso-8859-1">
|
|
<title>ATC: fix atc command</title>
|
|
<link href="doxygen.css" rel="stylesheet" type="text/css">
|
|
</head><body>
|
|
<!-- Generated by Doxygen 1.3.9.1 -->
|
|
<h1><a class="anchor" name="man_fix_atc">fix atc command</a></h1><h2><a class="anchor" name="syntax">
|
|
syntax</a></h2>
|
|
fix AtC transfer <type> <parameter_file><ul>
|
|
<li>type<br>
|
|
= thermal : thermal coupling with fields: temperature <br>
|
|
= two_temperature : electron-phonon coupling with field: temperature and electron_temperature <br>
|
|
= hardy : Hardy on-the-fly post-processing (see "related" section for possible fields) <br>
|
|
</li><li>parameter_file = name of the file with material parameters. <br>
|
|
note: hardy does not require a parameter file </li></ul>
|
|
<h2><a class="anchor" name="examples">
|
|
examples</a></h2>
|
|
<code> fix_modify AtC transfer thermal Ar_thermal.dat </code> <br>
|
|
<code> fix_modify AtC transfer hardy </code> <h2><a class="anchor" name="description">
|
|
description</a></h2>
|
|
This fix is the beginning to creating a coupled FE/MD simulation and/or an on-the-fly estimation of continuum fields. The coupled versions of this fix do Verlet integration and the Hardy/post-processing does not. After instantiating this fix, several other fix_modify commands will be needed to set up the problem, e.g. define the finite element mesh and prescribe initial and boundary conditions.<p>
|
|
The following coupling example is typical, but non-exhaustive:<br>
|
|
<p>
|
|
<code> # ... commands to create and initialize the MD system <br>
|
|
</code><p>
|
|
<code> # initial fix to designate coupling type and group to apply it to <br>
|
|
# tag group physics material_file <br>
|
|
fix AtC internal atc thermal Ar_thermal.mat<br>
|
|
<br>
|
|
# create a uniform 12 x 2 x 2 mesh that covers region contain the group <br>
|
|
# nx ny nz region periodicity <br>
|
|
fix_modify AtC fem create mesh 12 2 2 mdRegion f p p<br>
|
|
<br>
|
|
# specify the control method for the type of coupling <br>
|
|
# physics control_type <br>
|
|
fix_modify AtC transfer thermal control flux <br>
|
|
<br>
|
|
# specify the initial values for the empirical field "temperature" <br>
|
|
# field node_group value <br>
|
|
fix_modify AtC transfer initial temperature all 30.<br>
|
|
<br>
|
|
# create an output stream for nodal fields <br>
|
|
# filename output_frequency <br>
|
|
fix_modify AtC transfer output atc_fe_output 100<br>
|
|
<br>
|
|
</code><p>
|
|
<code> run 1000 <br>
|
|
</code><p>
|
|
likewise for this post-processing example: <br>
|
|
<p>
|
|
<code> # ... commands to create and initialize the MD system <br>
|
|
</code><p>
|
|
<code> # initial fix to designate post-processing and the group to apply it to <br>
|
|
# no material file is allowed nor required <br>
|
|
fix AtC internal atc hardy <br>
|
|
<br>
|
|
# create a uniform 1 x 1 x 1 mesh that covers region contain the group <br>
|
|
# with periodicity this effectively creats a system average <br>
|
|
fix_modify AtC fem create mesh 1 1 1 box p p p <br>
|
|
<br>
|
|
# change from default lagrangian map to eulerian <br>
|
|
# refreshed every 100 steps <br>
|
|
fix_modify AtC atom_element_map eulerian 100 <br>
|
|
<br>
|
|
# start with no field defined <br>
|
|
fix_modify AtC transfer fields none <br>
|
|
<br>
|
|
# add mass density, potential energy density, stress and temperature <br>
|
|
fix_modify AtC transfer fields add density energy stress temperature <br>
|
|
<br>
|
|
# create an output stream for nodal fields <br>
|
|
# filename output_frequency <br>
|
|
fix_modify AtC transfer output nvtFE 100 text <br>
|
|
</code><p>
|
|
<code> run 1000 <br>
|
|
</code><p>
|
|
Note coupling and post-processing can be combined in the same simulations using separate fixes. <br>
|
|
For detailed exposition of the theory and algorithms please see:<br>
|
|
<ul>
|
|
<li>Wagner, GJ; Jones, RE; Templeton, JA; Parks, MA. <em> An atomistic-to-continuum coupling method for heat transfer in solids. </em> Special Issue of Computer Methods and Applied Mechanics (2008) 197:3351.</li><li>Zimmerman, JA; Webb, EB; Hoyt, JJ;. Jones, RE; Klein, PA; Bammann, DJ, <em> Calculation of stress in atomistic simulation </em> Special Issue of Modelling and Simulation in Materials Science and Engineering (2004), 12:S319</li></ul>
|
|
<p>
|
|
Please refer to the standard finite element (FE) texts, e.g. T.J.R Hughes <em> The finite element method </em>, Dover 2003, for the basics of FE simulation.<h2><a class="anchor" name="restrictions">
|
|
restrictions</a></h2>
|
|
Thermal and two_temperature (coupling) types use a Verlet time-integration algorithm. The hardy type does not contain its own time-integrator and must be used with a separate fix that does contain one, e.g. nve, nvt, etc.<p>
|
|
Currently,<ul>
|
|
<li>the coupling is restricted to thermal physics</li><li>the FE computations are done in serial on each processor.</li></ul>
|
|
<h2><a class="anchor" name="related">
|
|
related</a></h2>
|
|
fix_modify commands for setup: <br>
|
|
<ul>
|
|
<li><a class="el" href="man_fem_mesh.html">fix_modify AtC fem create mesh</a></li><li><a class="el" href="man_mesh_nodeset.html">fix_modify AtC mesh create_nodeset</a></li><li><a class="el" href="man_mesh_faceset.html">fix_modify AtC mesh create_faceset</a></li><li><a class="el" href="man_mesh_elemset.html">fix_modify AtC mesh create_elementset</a></li><li><a class="el" href="man_transfer_internal.html">fix_modify AtC transfer internal</a></li><li><a class="el" href="man_transfer_boundary.html">fix_modify AtC transfer boundary</a></li><li><a class="el" href="man_internal_quadrature.html">fix_modify AtC transfer internal_quadrature</a></li><li><a class="el" href="man_time_integration.html">fix_modify AtC transfer pmfc</a></li><li><a class="el" href="man_electron_integration.html">fix_modify AtC extrinsic electron_integration</a></li></ul>
|
|
<p>
|
|
fix_modify commands for boundary and initial conditions:<br>
|
|
<ul>
|
|
<li><a class="el" href="man_initial.html">fix_modify AtC transfer initial</a></li><li><a class="el" href="man_fix_nodes.html">fix_modify AtC transfer fix</a></li><li><a class="el" href="man_unfix_nodes.html">fix_modify AtC transfer unfix</a></li><li><a class="el" href="man_fix_flux.html">fix_modify AtC transfer fix_flux</a></li><li><a class="el" href="man_unfix_flux.html">fix_modify AtC transfer unfix_flux</a></li><li><a class="el" href="man_source.html">fix_modify AtC transfer source</a></li><li><a class="el" href="man_remove_source.html">fix_modify AtC transfer remove_source</a></li></ul>
|
|
<p>
|
|
fix_modify commands for control and filtering: <br>
|
|
<ul>
|
|
<li><a class="el" href="man_thermal_control.html">fix_modify AtC transfer thermal control</a></li><li><a class="el" href="man_time_filter.html">fix_modify AtC transfer filter</a></li><li><a class="el" href="man_filter_scale.html">fix_modify AtC transfer filter scale</a></li><li><a class="el" href="man_equilibrium_start.html">fix_modify AtC transfer equilibrium_start</a></li><li><a class="el" href="man_extrinsic_exchange.html">fix_modify AtC extrinsic exchange</a></li></ul>
|
|
<p>
|
|
fix_modify commands for output: <br>
|
|
<ul>
|
|
<li><a class="el" href="man_transfer_output.html">fix_modify AtC transfer output</a></li><li><a class="el" href="man_transfer_atomic_output.html">fix_modify AtC transfer atomic_output</a></li><li><a class="el" href="man_mesh_output.html">fix_modify AtC mesh output</a></li><li><a class="el" href="man_write_restart.html">fix_modify AtC transfer write_restart</a></li><li><a class="el" href="man_read_restart.html">fix_modify AtC transfer read_restart</a></li></ul>
|
|
<p>
|
|
fix_modify commands for post-processing: <br>
|
|
<ul>
|
|
<li><a class="el" href="man_hardy_fields.html">fix_modify AtC transfer fields</a></li><li><a class="el" href="man_hardy_gradients.html">fix_modify AtC transfer gradients</a></li><li><a class="el" href="man_hardy_rates.html">fix_modify AtC transfer rates</a></li><li><a class="el" href="man_hardy_computes.html">fix_modify AtC transfer computes</a></li><li><a class="el" href="man_hardy_set.html">fix_modify AtC set</a></li><li><a class="el" href="man_hardy_on_the_fly.html">fix_modify AtC transfer on_the_fly</a></li><li><a class="el" href="man_boundary_integral.html">fix_modify AtC boundary_integral</a></li><li><a class="el" href="man_contour_integral.html">fix_modify AtC contour_integral</a></li></ul>
|
|
<p>
|
|
miscellaneous fix_modify commands: <br>
|
|
<ul>
|
|
<li><a class="el" href="man_atom_element_map.html">fix_modify AtC transfer atom_element_map</a></li><li><a class="el" href="man_neighbor_reset_frequency.html">fix_modify AtC transfer neighbor_reset_frequency</a></li></ul>
|
|
<p>
|
|
Note: a set of example input files with the attendant material files are included with this package <h2><a class="anchor" name="default">
|
|
default</a></h2>
|
|
none <hr size="1"><address style="align: right;"><small>Generated on Mon Aug 17 09:35:16 2009 for ATC by
|
|
<a href="http://www.doxygen.org/index.html">
|
|
<img src="doxygen.png" alt="doxygen" align="middle" border="0"></a> 1.3.9.1 </small></address>
|
|
</body>
|
|
</html>
|