lammps/doc/USER/atc/man_fix_atc.html

255 lines
15 KiB
HTML

<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<title>ATC: fix atc command</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<link href="doxygen.css" rel="stylesheet" type="text/css"/>
</head>
<body>
<!-- Generated by Doxygen 1.6.1 -->
<div class="navigation" id="top">
<div class="tabs">
<ul>
<li><a href="index.html"><span>Main&nbsp;Page</span></a></li>
<li class="current"><a href="pages.html"><span>Related&nbsp;Pages</span></a></li>
<li><a href="namespaces.html"><span>Namespaces</span></a></li>
<li><a href="annotated.html"><span>Classes</span></a></li>
<li><a href="files.html"><span>Files</span></a></li>
<li><a href="dirs.html"><span>Directories</span></a></li>
</ul>
</div>
</div>
<div class="contents">
<h1><a class="anchor" id="man_fix_atc">fix atc command </a></h1><h2><a class="anchor" id="syntax">
syntax</a></h2>
<p>fix &lt;fixID&gt; &lt;group&gt; atc &lt;type&gt; &lt;parameter_file&gt;</p>
<ul>
<li>fixID = name of fix</li>
<li>group = name of group fix is to be applied</li>
<li>type<br/>
= thermal : thermal coupling with fields: temperature <br/>
= two_temperature : electron-phonon coupling with field: temperature and electron_temperature <br/>
= hardy : on-the-fly post-processing using kernel localization functions (see "related" section for possible fields) <br/>
= field : on-the-fly post-processing using mesh-based localization functions (see "related" section for possible fields) <br/>
</li>
<li>parameter_file = name of the file with material parameters. <br/>
note: Neither hardy nor field requires a parameter file </li>
</ul>
<h2><a class="anchor" id="examples">
examples</a></h2>
<p><code> fix AtC internal atc thermal Ar_thermal.dat </code> <br/>
<code> fix AtC internal atc two_temperature Ar_ttm.mat </code> <br/>
<code> fix AtC internal atc hardy </code> <br/>
<code> fix AtC internal atc field </code> <br/>
</p>
<h2><a class="anchor" id="description">
description</a></h2>
<p>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 /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>
<p>The following coupling example is typical, but non-exhaustive:<br/>
</p>
<p><code> # ... commands to create and initialize the MD system <br/>
</code></p>
<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 mesh create 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 thermal control flux <br/>
<br/>
# specify the initial values for the empirical field "temperature" <br/>
# field node_group value <br/>
fix_modify AtC initial temperature all 30.<br/>
<br/>
# create an output stream for nodal fields <br/>
# filename output_frequency <br/>
fix_modify AtC output atc_fe_output 100<br/>
<br/>
</code></p>
<p><code> run 1000 <br/>
</code></p>
<p>likewise for this post-processing example: <br/>
</p>
<p><code> # ... commands to create and initialize the MD system <br/>
</code></p>
<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/>
# for hardy fix, specific kernel function (function type and range) to # be used as a localization function <br/>
fix AtC kernel quartic_sphere 10.0 <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 mesh create 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/>
# add mass density, potential energy density, stress and temperature <br/>
fix_modify AtC fields add density energy stress temperature <br/>
<br/>
# create an output stream for nodal fields <br/>
# filename output_frequency <br/>
fix_modify AtC output nvtFE 100 text <br/>
</code></p>
<p><code> run 1000 <br/>
</code></p>
<p>the mesh's linear interpolation functions can be used as the localization function <br/>
by using the field option: <br/>
</p>
<p><code> fix AtC internal atc field <br/>
<br/>
fix_modify AtC mesh create 1 1 1 box p p p <br/>
<br/>
... <br/>
<br/>
</code></p>
<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/>
</p>
<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. <br/>
</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. <br/>
</li>
<li>Zimmerman, JA; Jones, RE; Templeton, JA, <em> A material frame approach for evaluating continuum variables in atomistic simulations. </em> Journal of Computational Physics (2010), 229:2364. <br/>
</li>
<li>Templeton, JA; Jones, RE; Wagner, GJ, <em> Application of a field-based method to spatially varying thermal transport problems in molecular dynamics. </em> Modelling and Simulation in Materials Science and Engineering (2010), 18:085007. <br/>
</li>
<li>Jones, RE; Templeton, JA; Wagner, GJ; Olmsted, D; Modine, JA, <em> Electron transport enhanced molecular dynamics for metals and semi-metals. </em> International Journal for Numerical Methods in Engineering (2010), 83:940. <br/>
</li>
<li>Templeton, JA; Jones, RE; Lee, JW; Zimmerman, JA; Wong, BM, <em> A long-range electric field solver for molecular dynamics based on atomistic-to-continuum modeling. </em> Journal of Chemical Theory and Computation (2011), 7:1736. <br/>
</li>
<li>Mandadapu, KK; Templeton, JA; Lee, JW, <em> Polarization as a field variable from molecular dynamics simulations. </em> Journal of Chemical Physics (2013), 139:054115. <br/>
</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.</p>
<h2><a class="anchor" id="restrictions">
restrictions</a></h2>
<p>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>
<p>Currently,</p>
<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" id="related">
related</a></h2>
<p>fix_modify commands for setup: <br/>
</p>
<ul>
<li><a class="el" href="man_mesh_create.html">fix_modify AtC mesh create</a></li>
<li><a class="el" href="man_mesh_quadrature.html">fix_modify AtC mesh quadrature</a></li>
<li><a class="el" href="man_mesh_read.html">fix_modify AtC mesh read</a></li>
<li><a class="el" href="man_mesh_write.html">fix_modify AtC mesh write</a></li>
<li><a class="el" href="man_mesh_create_nodeset.html">fix_modify AtC mesh create_nodeset</a></li>
<li><a class="el" href="man_mesh_add_to_nodeset.html">fix_modify AtC mesh add_to_nodeset</a></li>
<li><a class="el" href="man_mesh_create_faceset_box.html">fix_modify AtC mesh create_faceset box</a></li>
<li><a class="el" href="man_mesh_create_faceset_plane.html">fix_modify AtC mesh create_faceset plane</a></li>
<li><a class="el" href="man_mesh_create_elementset.html">fix_modify AtC mesh create_elementset</a></li>
<li><a class="el" href="man_mesh_delete_elements.html">fix_modify AtC mesh delete_elements</a></li>
<li><a class="el" href="man_mesh_nodeset_to_elementset.html">fix_modify AtC mesh nodeset_to_elementset</a></li>
<li><a class="el" href="man_boundary.html">fix_modify AtC boundary</a></li>
<li><a class="el" href="man_internal_quadrature.html">fix_modify AtC internal_quadrature</a></li>
<li><a class="el" href="man_thermal_time_integration.html">fix_modify AtC time_integration (thermal)</a></li>
<li><a class="el" href="man_momentum_time_integration.html">fix_modify AtC time_integration (momentum)</a></li>
<li><a class="el" href="man_electron_integration.html">fix_modify AtC extrinsic electron_integration</a></li>
<li><a class="el" href="man_internal_element_set.html">fix_modify AtC internal_element_set</a></li>
<li><a class="el" href="man_decomposition.html">fix_modify AtC decomposition</a></li>
</ul>
<p>fix_modify commands for boundary and initial conditions:<br/>
</p>
<ul>
<li><a class="el" href="man_initial.html">fix_modify AtC initial</a></li>
<li><a class="el" href="man_fix_nodes.html">fix_modify AtC fix</a></li>
<li><a class="el" href="man_unfix_nodes.html">fix_modify AtC unfix</a></li>
<li><a class="el" href="man_fix_flux.html">fix_modify AtC fix_flux</a></li>
<li><a class="el" href="man_unfix_flux.html">fix_modify AtC unfix_flux</a></li>
<li><a class="el" href="man_source.html">fix_modify AtC source</a></li>
<li><a class="el" href="man_remove_source.html">fix_modify AtC remove_source</a></li>
</ul>
<p>fix_modify commands for control and filtering: <br/>
</p>
<ul>
<li><a class="el" href="man_control.html">fix_modify AtC control</a></li>
<li><a class="el" href="man_control_thermal.html">fix_modify AtC control thermal</a></li>
<li><a class="el" href="man_control_thermal_correction_max_iterations.html">fix_modify AtC control thermal correction_max_iterations</a></li>
<li><a class="el" href="man_control_momentum.html">fix_modify AtC control momentum</a></li>
<li><a class="el" href="man_localized_lambda.html">fix_modify AtC control localized_lambda</a></li>
<li><a class="el" href="man_lumped_lambda_solve.html">fix_modify AtC control lumped_lambda_solve</a></li>
<li><a class="el" href="man_mask_direction.html">fix_modify AtC control mask_direction</a></li>
<li><a class="el" href="man_time_filter.html">fix_modify AtC filter</a></li>
<li><a class="el" href="man_filter_scale.html">fix_modify AtC filter scale</a></li>
<li><a class="el" href="man_filter_type.html">fix_modify AtC filter type</a></li>
<li><a class="el" href="man_equilibrium_start.html">fix_modify AtC equilibrium_start</a></li>
<li><a class="el" href="man_extrinsic_exchange.html">fix_modify AtC extrinsic exchange</a></li>
<li><a class="el" href="man_poisson_solver.html">fix_modify AtC poisson_solver</a></li>
</ul>
<p>fix_modify commands for output: <br/>
</p>
<ul>
<li><a class="el" href="man_output.html">fix_modify AtC output</a></li>
<li><a class="el" href="man_output_nodeset.html">fix_modify AtC output nodeset</a></li>
<li><a class="el" href="man_output_elementset.html">fix_modify AtC output elementset</a></li>
<li><a class="el" href="man_boundary_integral.html">fix_modify AtC output boundary_integral</a></li>
<li><a class="el" href="man_contour_integral.html">fix_modify AtC output contour_integral</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 write_restart</a></li>
<li><a class="el" href="man_read_restart.html">fix_modify AtC read_restart</a></li>
</ul>
<p>fix_modify commands for post-processing: <br/>
</p>
<ul>
<li><a class="el" href="man_hardy_kernel.html">fix_modify AtC kernel</a></li>
<li><a class="el" href="man_hardy_fields.html">fix_modify AtC fields</a></li>
<li><a class="el" href="man_hardy_gradients.html">fix_modify AtC gradients</a></li>
<li><a class="el" href="man_hardy_rates.html">fix_modify AtC rates</a></li>
<li><a class="el" href="man_hardy_computes.html">fix_modify AtC computes</a></li>
<li><a class="el" href="man_hardy_on_the_fly.html">fix_modify AtC on_the_fly</a></li>
<li><a class="el" href="man_pair_interactions.html">fix_modify AtC pair_interactions/bond_interactions</a></li>
<li><a class="el" href="man_sample_frequency.html">fix_modify AtC sample_frequency</a></li>
<li><a class="el" href="man_set.html">fix_modify AtC set</a></li>
</ul>
<p>miscellaneous fix_modify commands: <br/>
</p>
<ul>
<li><a class="el" href="man_atom_element_map.html">fix_modify AtC atom_element_map</a></li>
<li><a class="el" href="man_atom_weight.html">fix_modify AtC atom_weight</a></li>
<li><a class="el" href="man_write_atom_weights.html">fix_modify AtC write_atom_weights</a></li>
<li><a class="el" href="man_reset_time.html">fix_modify AtC reset_time</a></li>
<li><a class="el" href="man_reset_atomic_reference_positions.html">fix_modify AtC reset_atomic_reference_positions</a></li>
<li><a class="el" href="man_fe_md_boundary.html">fix_modify AtC fe_md_boundary</a></li>
<li><a class="el" href="man_boundary_faceset.html">fix_modify AtC boundary_faceset</a></li>
<li><a class="el" href="man_consistent_fe_initialization.html">fix_modify AtC consistent_fe_initialization</a></li>
<li><a class="el" href="man_mass_matrix.html">fix_modify AtC mass_matrix</a></li>
<li><a class="el" href="man_material.html">fix_modify AtC material</a></li>
<li><a class="el" href="man_atomic_charge.html">fix_modify AtC atomic_charge</a></li>
<li><a class="el" href="man_source_integration.html">fix_modify AtC source_integration</a></li>
<li><a class="el" href="man_temperature_definition.html">fix_modify AtC temperature_definition</a></li>
<li><a class="el" href="man_track_displacement.html">fix_modify AtC track_displacement</a></li>
<li><a class="el" href="man_boundary_dynamics.html">fix_modify AtC boundary_dynamics</a></li>
<li><a class="el" href="man_add_species.html">fix_modify AtC add_species</a></li>
<li><a class="el" href="man_add_molecule.html">fix_modify AtC add_molecule</a></li>
<li><a class="el" href="man_remove_species.html">fix_modify AtC remove_species</a></li>
<li><a class="el" href="man_remove_molecule.html">fix_modify AtC remove_molecule</a></li>
</ul>
<p>Note: a set of example input files with the attendant material files are included with this package </p>
<h2><a class="anchor" id="default">
default</a></h2>
<p>none </p>
</div>
<hr size="1"/><address style="text-align: right;"><small>Generated on 21 Aug 2013 for ATC by&nbsp;
<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.6.1 </small></address>
</body>
</html>