File: workedExample.tex

package info (click to toggle)
xmds-doc 0~svn.1884-3.1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd, wheezy
  • size: 8,336 kB
  • ctags: 192
  • sloc: makefile: 135; python: 55
file content (715 lines) | stat: -rw-r--r-- 30,556 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
% $Id: workedExample.tex 1339 2007-05-09 05:07:30Z joehope $

% Copyright (C) 2000-2007
%
% Code contributed by Greg Collecutt, Joseph Hope and the xmds-devel team
%
% This file is part of xmds.
%
% This program is free software; you can redistribute it and/or
% modify it under the terms of the GNU General Public License
% as published by the Free Software Foundation; either version 2
% of the License, or (at your option) any later version.
%
% This program is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with this program; if not, write to the Free Software 
% Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.

\chapter{Worked example: nlse.xmds} 
\label{chap:exampleNlse}

XML is very much like HTML, so readers with any experience in writing
HTML will find XML fairly straight forward. However, even if a user
has not had any experience in HTML they will pick up the XML with the
help of some example simulation files. A comprehensive guide to
general XML syntax is available in~\cite{web:XML_syntax}. For now, a
good introduction will be to step through the ``\ttt{nlse.xmds}''
example script, then we can work through some more advanced examples
to elaborate on some of the options available. These examples scripts can 
be found in the \ttt{/examples} directory, along with many others. They may 
also be downloaded from the web site
\htmladdnormallink{http://www.xmds.org/examples.html}{http://www.xmds.org/examples.html}.

\begin{xmdsCode}
<?xml version="1.0"?>
<!--Example simulation: Non-Linear Schroedinger Equation-->

<simulation>

  <name>nlse</name>
  <prop_dim>z</prop_dim>
  <error_check>yes</error_check>
  <stochastic>no</stochastic>

  <globals>
    <![CDATA[
      const double energy = 4;
      const double vel = 0.0;
      const double hwhm = 1.0;
    ]]>
  </globals>

  <field>
    <name>main</name>
    <dimensions>  t   </dimensions>
    <lattice>    100  </lattice>
    <domains>  (-5,5) </domains>

    <samples>1</samples>
    <vector>
      <name>main</name>
      <type>complex</type>
      <components>phi</components>
      <fourier_space>no</fourier_space>
      <![CDATA[
        const double w0 = hwhm*sqrt(2/log(2));
        const double amp  = sqrt(energy/w0/sqrt(M_PI/2));

        phi = pcomplex(amp*exp(-t*t/w0/w0),vel*t);
      ]]>
    </vector>
    <vector>
      <name>vc1</name>
      <type>double</type>
      <components>damping</components>
      <fourier_space>no</fourier_space>
      <![CDATA[
        damping=1.0*(1-exp(-pow(t*t/4/4,10)));
      ]]>
    </vector>

  </field>

  <sequence>
    <integrate>

      <algorithm>RK4IP</algorithm>
      <interval>20</interval>
      <lattice>1000</lattice>
      <samples>50</samples>
      <k_operators>
        <constant>yes</constant>
        <operator_names>L</operator_names>
        <![CDATA[
          L = rcomplex(0,-kt*kt/2);
      ]]>
      </k_operators>
      <vectors>main vc1</vectors>
      <![CDATA[
        dphi_dz =  L[phi] + i*~phi*phi*phi - phi*damping;
      ]]>
    </integrate>
  </sequence>

  <output>
    <filename>nlse.xsil</filename>
    <group>
      <sampling>
        <fourier_space> no </fourier_space>
        <lattice>       50 </lattice>
        <moments>pow_dens</moments>
        <![CDATA[
          pow_dens=~phi*phi;
        ]]>
      </sampling>
    </group>
  </output>

</simulation>
\end{xmdsCode}

The above XML file describes to \xmds how to write a program that
solves the Nonlinear Schr\"{o}dinger Equation in one dimension;
\eqn{eq:nlse1d}.
\begin{equation}
\frac{\partial \phi}{\partial \xi } = i\left[\frac{1}{2} \frac{\partial ^{2}
 \phi}{\partial \tau ^{2}} + i \Gamma (\tau) \phi + |\phi|^{2} \phi
 \right].
\label{eq:nlse1d}
\end{equation}
Here $\phi$ is the single component complex field (\ttt{phi} in the
input file), $\xi$ is the propagation dimension (\ttt{z}), $\tau$ is
the local time coordinate (\ttt{t}), and $\Gamma$ is a damping field
(\ttt{damping}) applied to absorb scattered radiation at the domain
boundaries.  \fig{fig:nlse} shows what the output of this
simulation should look like.
\begin{figure}[ht]
\centerline{\includegraphics[width=\figwidth]{figures/nlse}}
\caption{Results for nlse.xmds}
\label{fig:nlse}
\end{figure}

At the head of the file:
\begin{xmdsCode}
<?xml version="1.0"?>
<!--Example simulation: Non-Linear Schroedinger Equation-->
\end{xmdsCode}
Here the first line defines the file as an XML file. \xmds
requires correct XML syntax, but otherwise does not need to validate
the file with a DTD (Document Type Definition). The second line is
simply a comment line. All comments begin with a ``\ttt{<!--}'' and
end with a ``\ttt{-->}''. Any characters are allowed in-between save
for the comment closing sequence ``\ttt{-->}''.

The remainder of the file is covered in detail in
Sections~\ref{sec:simulationElement}--\ref{sec:outputElement}, while
\Sec{sec:syntaxSummary} contains a tabulated summary of the syntax.

\section{The simulation element} 
\label{sec:simulationElement}

\begin{xmdsCode}
<simulation>

  <name>nlse</name>
  <prop_dim>z</prop_dim>
  <error_check>yes</error_check>
  <stochastic>no</stochastic>
  
  <!-- More xmds code -->

</simulation>
\end{xmdsCode}

The root XML element of the input file must be a \xmdsTag{simulation}
element, which contains the high level description of the problem.

Within the \xmdsTag{simulation} element \xmds will first look for
a simulation \xmdsTag{name}. This information is optional, the default
being the file name minus its last extension. It is used as the base
for the output \ttt{.cc} file, and is the name of the executable that
\xmds will make.

The \xmdsTag{prop\_dim} assignment defines the name of the main
propagation dimension. It is compulsory, as it is expected that
whichever symbol is used here will also appear in the equations the
user supplies.

The \xmdsTag{error\_check} assignment is optional, the default being
``\ttt{yes}''. This option tells \xmds whether to write the code
so as to repeat the simulation with half the step sizes in the
integration segments, and subsequently compare the two sets of
results. It is a good idea to leave this option on until confident
that the errors are acceptable for the simulation in question. If set
to ``\ttt{yes}'' the output executable will compare all output
moments on a point by point basis between the full and half-step runs.
It will then write out to screen the maximum value of the discrepancy,
and the half-step results will be used for the final output.

The \xmdsTag{stochastic} assignment is also optional, the default being
``\ttt{yes}''. It tells \xmds whether the simulation uses Gaussian
noise terms or not. If \xmdsTag{stochastic} is set to ``\ttt{yes}''
then three further assignments, \xmdsTag{paths}, \xmdsTag{seed}, and
\xmdsTag{noises} become compulsory. The \ttt{kubo.xmds} and
\ttt{fibre.xmds} examples in \Chap{chap:moreExamples} cover stochastic
simulations in more detail.

The child elements within the \xmdsTag{simulation} element that \xmds
then looks for are \xmdsTag{globals}, \xmdsTag{field},
\xmdsTag{sequence}, and \xmdsTag{output}.

\section{The globals element} 
\label{sec:globalsElement}

\begin{xmdsCode}
<globals>
  <![CDATA[
    const double energy = 4;
    const double vel = 0.0;
    const double hwhm = 1.0;
  ]]>
</globals>
\end{xmdsCode}
This element is optional, and is used to define any numerical
constants that are useful to have globally available to all sections
of code. The \ttt{<![CDATA[ ... ]]>} container tells the XML parser
to copy its contents without attempting to parse them. The contents:
\begin{xmdsCode}
const double energy = 4;
const double vel = 0.0;
const double hwhm = 1.0;
\end{xmdsCode}
are copied directly into the output C code and therefore must have
correct C syntax. See \Sec{sec:c-coding} for more on C syntax.

\section{The field element}
\label{sec:fieldElement}

\begin{xmdsCode}
<field>
  <name>main</name>
  <dimensions>  t   </dimensions>
  <lattice>    100  </lattice>
  <domains>  (-5,5) </domains>

  <samples>1</samples>
    
  <!-- Possibly more xmds code -->

</field>
\end{xmdsCode}

This element is compulsory as it is central to the problem definition:
it specifies the geometry of the field. There are five assignments
that \xmds looks for in this element: \xmdsTag{name},
\xmdsTag{dimensions}, \xmdsTag{lattice}, \xmdsTag{domains}, and
\xmdsTag{samples}.

The \xmdsTag{name} assignment provides a name for the field. It is not
compulsory and will default to ``\ttt{main}''. It is envisaged that
future versions \xmds may allow sub fields within fields, but at
present only one \xmdsTag{field} element is permitted.

The \xmdsTag{dimensions} assignment lists the names of the transverse
field dimensions. From this point on the number of transverse
dimensions is set to the number of names found in this assignment. If
this assignment is empty or absent then the field is assumed to be
without transverse dimensions (as in the \ttt{kubo.xmds} example),
and \xmds will not look for either a \xmdsTag{lattice} assignment
or a \xmdsTag{domains} assignment. The dimension names must be valid
names as far as the C language is concerned--see \Sec{sec:c-coding}.

The \xmdsTag{lattice} assignment is compulsory if there were one or
more transverse dimension names found in the \xmdsTag{dimensions}
assignment, and \xmds will look for the same number of positive
integers here as there were transverse dimensions. Each of these
integers defines the number of points or steps to use in the
corresponding transverse dimension, and so \xmds expects to find
integers with a value of 2 or more in each case.

The \xmdsTag{domains} assignment is similar, except that it defines the
domain range for each dimension. This is done by entering a bracketed
pair of real numbers for each dimension, and \xmds expects the
range $b-a$ for any pair (a,b) to be greater than $10^{-100}$. Each
domain range is divided according to \eqn{eq:domRange}, in
which $n$ is the number of lattice points. Note that the points stop
one spacing short of the upper value---this is due to the forced
periodic boundary conditions.
\begin{equation}
x_{j} = \frac{(n-j+1)a + (j-1)b}{n}, \quad j=1,\dots,n.
\label{eq:domRange}
\end{equation}

Finally the \xmdsTag{samples} assignment tells \xmds which moment
groups (\Sec{sec:outputElement}) to sample immediately after
all the field's vectors have been initialised. A sequence of
``\ttt{0}'' or ``\ttt{1}'' entries is expected, as many as there are
moment groups defined in the output (\xmds checks this once all
moment groups have been processed). The correspondence is
sequential. A ``\ttt{1}'' means sample the moment group in question,
a ``\ttt{0}'' means do not.

Next \xmds will look for the field's vectors.

\section{The vector element} 
\label{sec:vectorElement}

\begin{xmdsCode}
<vector>
  <name>main</name>
  <type>complex</type>
  <components>phi</components>
  <fourier_space>no</fourier_space>
  <![CDATA[
    const double w0 = hwhm*sqrt(2/log(2));
    const double amp  = sqrt(energy/w0/sqrt(M_PI/2));

    phi = pcomplex(amp*exp(-t*t/w0/w0),vel*t);
  ]]>
</vector>
\end{xmdsCode}

The \xmdsTag{name} assignment is compulsory here, as it is quite common
to employ more than one vector. The field must have at least one
vector called ``\ttt{main}'', as it is this vector that \xmds
forward evolves in the main propagation dimension.

The \xmdsTag{type} assignment is optional. It will default to
``\ttt{complex}'' unless specifically stated as ``\ttt{double}'',
which is the standard data type in C for floating point
variables. Note that vectors of type ``\ttt{double}'' cannot be
Fourier transformed without some arbitrary definition of how this is
done. Thus \xmds requires that ``\ttt{double}'' vectors always
remain in the Fourier space in which they were initialised. A request
to access them in some other Fourier space will result in an error.

The \xmdsTag{components} assignment is compulsory. It identifies by
name the components of the vector. Here there is only one component,
but more may be defined simply by separating them with spaces, for
example \xmdsTag{components}\ttt{phi} \ttt{theta}\xmdsTag{/components}. The component names
must be valid as far as the C language is concerned -- see
\Sec{sec:c-coding}.

Finally the initialisation of the vector must be defined. This can be
done in two ways. Either by means of a \xmdsTag{filename} assignment,
or by supplying C code. Here we are doing the latter. In the absence
of \xmdsTag{filename} assignment \xmds looks for a
\xmdsTag{fourier\_space} assignment followed by the initialisation code.

The \xmdsTag{fourier\_space} assignment must contain as many
``\ttt{yes}'' or ``\ttt{no}'' entries as there are transverse
dimensions in the field, which in this case is only one. A
``\ttt{yes}'' entry means the corresponding dimension is in Fourier
space when initialised, and vice-versa.

When initialising from code, the code is included as the content of
the \xmdsTag{vector} element. There must be an equation for each
component, and these may include functions of \xmdsTag{global}
variables and of the transverse dimensions. Again see
\Sec{sec:c-coding} for more information. Two points to remember:
firstly if the field components are complex variables they may be
initialised using the complex functions in \tab{tab:functions}, and
secondly every field component should be explicitly initialised --
even if it is only being initialised to zero.

Finally, when initialising from file, \xmds expects that the file
is simply a sequence of numbers in ASCII format, with as many numbers
as the product of lattice points and vector components. In usual C
convention, the component index varies most rapidly, then the right
most transverse dimension lattice and so on. The left most transverse
dimension lattice dimension varies most slowly. If the vector is
complex then each data point is read in as a sequential real and
imaginary pair. All vectors initialised from file are assumed to be
entirely in normal space.

\section{The sequence element} 
\label{sec:sequenceElement}

\begin{xmdsCode}
<sequence>

  <integrate>
      
    <!-- More xmds code -->

  </integrate>
    
  <!-- More xmds code -->

</sequence>
\end{xmdsCode}

Next comes the \xmdsTag{sequence} element, which can be used in two
different contexts. The first context (as used here) is the ``top
level'' \xmdsTag{sequence} element, and is simply for defining the
sequence of segments that happen to the field after initialisation. In
a stochastic simulation this sequence is repeated for each integration
path, but in this non-stochastic example the sequence is only stepped
through once. In this context \xmds expects to find any number of
\xmdsTag{integrate}, \xmdsTag{filter}, or \xmdsTag{sequence} elements --
which leads us to the second context.

In the second context \xmdsTag{sequence} elements may be used as child
elements within parent \xmdsTag{sequence} elements, in which case they
represent a sub-loop to be repeated one or more times. The number of
repetitions is defined by a simple \xmdsTag{cycles} assignment within,
which is expected to be the first element. In the absence of a
\xmdsTag{cycles} element the number of cycles defaults to one, but if
it is present it must come first. Any number of \xmdsTag{integrate},
\xmdsTag{filter}, or \xmdsTag{sequence} elements may then
follow. Obviously the ordering of segments within any sequence element
is important.

\section{The integrate element} 
\label{sec:integrateElement}

\begin{xmdsCode}
<integrate>

  <algorithm>RK4IP</algorithm>
  <interval>20</interval>
  <lattice>1000</lattice>
  <samples>50</samples>
  <k_operators>
    <constant>yes</constant>
    <operator_names>L</operator_names>
    <![CDATA[
          L = rcomplex(0,-kt*kt/2);
    ]]>
  </k_operators>
  <vectors>main vc1</vectors>
  <![CDATA[
    dphi_dz =  L[phi] + i*~phi*phi*phi - phi*damping;
  ]]>
</integrate>
\end{xmdsCode}

Here lies the heart of what \xmds is all about, and not
surprisingly it forms the most complex element. The \xmdsTag{algorithm}
assignment is optional, and will default to \ttt{SIEX} for stochastic
simulations and \ttt{RK4EX} for non-stochastic simulations.  The
assignments \xmdsTag{interval}, \xmdsTag{lattice}, and \xmdsTag{samples}, are
all compulsory, and respectively represent the integration range,
total number of steps, and number of samples for each output moment
group to take within these steps. Thus each integer in the
\xmdsTag{samples} assignment must be either zero or else a factor of
\xmdsTag{lattice}, and if not \xmds will exit with an error
message to that effect.

The rest of the \xmdsTag{integrate} element consists of ``writing
down'' the form of \eqn{eq:xmdsPdeEx} in a high level format.  Here
the algorithm is \ttt{RK4IP}, which, as explained earlier, is an
interaction picture algorithm. The dispersion term in \eqn{eq:nlse1d}
is represented here by the linear operator ``\ttt{L}'' acting on
field component ``\ttt{phi}''. The presence of the term
``\ttt{L[phi]}'' in the equations causes the field to be evolved in
Fourier space in accordance with the algorithm listed in
\Sec{sec:rk4ipMethod}. The $\vect{\mathcal{N}}$ operators are
everything left over once terms such as this are removed.

The \xmdsTag{k\_operators} element contains just what its name implies:
k-space operators. Within this element a \xmdsTag{constant} declaration
is used to indicate whether or not these operators depend on the
propagation dimension. This declaration is optional with the default
being ``\ttt{no}'', but this will produce slower code. If the
operators do not depend on the propagation dimension then you will get
faster code is you enter \xmdsTag{constant}\ttt{yes}\xmdsTag{/constant}. An
\xmdsTag{operator\_names} assignment is compulsory. One or more operator
names may be listed here, and they all must be explicitly defined in
the code. Refer to \eqn{eq:l_maps} for how to write the
$\mathcal{L}^{ij}$ operators as functions of the Fourier space
dimensions.  The operator names must be valid names as far as the C
language is concerned---see \Sec{sec:c-coding}.

If there are no non-zero linear operators $\mathcal{L}^{ij}$, then the
entire \xmdsTag{k\_operators} element may be omitted.

\centerline{\bf Extreme DANGER!!!!!!!}

Because of the way \xmds works internally, it is possible when
using an interaction picture routine to write the main equations in
such as way as that everything works fine, but the results will be
incorrect! The problem is that \xmds does not closely inspect the
main equations. What it does is it looks for operator[component]
combinations, registers their presence, and the replaces their text
with something else. In the explicit picture they are replaced with a
reference to an internally calculated vector ($\vect{p} =
\mathcal{F}^{-1}\left[\mathcal{L}\left(x^0,\vect{k_\bot}\right) \cdot 
\mathcal{F}\left[\vect{a}\right]\right]$),
but in the interaction picture they are replaced with ``0'' and the
component $a^i$ being acted on by the operator $\mathcal{L}$ is
evolved in Fourier space using $a^i =
e^{\frac{1}{2}h\mathcal{L}^i\left(x^0,\vect{k_\bot}\right)} a^i$. This
means that if you write
\begin{xmdsCode}
<![CDATA[
  dtheta_dz =  L1[theta] + i*theta*theta*theta - theta*damping;
  dphi_dz   =  L2[theta] + i*~phi *phi  *phi   - phi  *damping; // no !!
]]>
\end{xmdsCode}
the component \ttt{theta} will be evolved in Fourier space using the
sum of the operators \ttt{L1} and \ttt{L2}, and the evolution of
\ttt{phi} will not include \ttt{L2[theta]}. Further, writing
something like \ttt{3*L1[theta]} does not have the intended effect -
\ttt{theta} will still only be evolved with \ttt{L1}, not
\ttt{3*L1}. We realise this is a potential source of error, but it
enables the equation syntax to remain uniform for all algorithms. A
future version of \xmds will check the users equations more
thoroughly.  When using explicit picture algorithms none of these
problems exist. If no errors are reported then what you get will be
what your equations asked for.

Sometimes, it is useful to calculate a function without having to calculate it separately for every 
transverse position.  This can be done by including the code in a \xmdsTag{functions} element.
\begin{xmdsCode}
<functions>
  <![CDATA[
  	double f = sin(2.0*z);  // This is a function of the propagation dimension only.
  ]]>
</functions>
\end{xmdsCode}
Conversely, it is often desirable to reference a variable in the equations that is
a function of the spatial coordinates, but is otherwise constant.  The
variable \ttt{damping} is exactly one of these. Rather than
re-calculating it at every time step and lattice point, what we did
here was to define an additional field vector ``\ttt{vc1}'' in the
\xmdsTag{field} element that could store this variable, calculate it
once when it is initialised, and then reference it in the main
equations. Therefore, we must tell \xmds which vectors we wish to
access in the equations. This is done with the \xmdsTag{vectors}
assignment. The reason that we do not make
{\em all} vectors available everywhere by default is one of
efficiency. For example if a particular vector was initialised in
x-space then it would have to be transformed to k-space in order to
use it in the \xmdsTag{k\_operators} code. If we did not need to use it
there then the Fourier transform would have been a waste of CPU
time. Further, the \ttt{damping} variable used here would have to
have been declared as type \ttt{complex} in order to have been able
to be Fourier transformed to k-space to be available in the
\xmdsTag{k\_operators} element where it wasn't needed, causing a waste of
memory as well.

Another possible function required to define the equations of motion is some integral across the transverse dimensions.  These are specified with the \xmdsTag{moment\_group} element.  We do not have one of these in our worked example, but as an example we might need the total number of atoms:
\begin{xmdsCode}
   <moment_group>
       <moments>number</moments>
       <integrate_dimension>yes</integrate_dimension>
	   <![CDATA[
		number += ~phi*phi;   
	]]>
      </moment_group>    
\end{xmdsCode}
The child elements here name the moments to be defined, and then describe which, if any, of the transverse dimensions are to be integrated.  The \ttt{CDATA} block defines the moments themselves.  The "number" moment in this example will be integrated over the only transverse dimension.  \xmdsTag{moment\_group} and \xmdsTag{functions} elements may be used any number of times in any order.  The propagation of the integration equations themselves is performed where the \xmdsTag{vectors} tag is placed.

Something else to keep in mind is that when a smooth function is
required it is often tempting to use transcendental functions (such as
\ttt{sqrt}, \ttt{exp}, \ttt{sin}, \ttt{cos}, \ttt{log}, ...). However,
these functions are computationally expensive, particularly on
processors that were not designed with such functions in mind. If they
are only used once to pre-calculate a constant vector (as is done
here) then fine, but if you include them in your main equations then
expect a big reduction in performance. If a smooth function that
depends on the propagation dimension is required, then, unless it is
essential to the model, it is usually better use a low order
polynomial approximation.

\section{The filter element} 
\label{sec:filterElement}

The \xmdsTag{filter} element, though not covered in this example, is
relatively straightforward. Here is an example which retains a
Gaussian profile of the field centered on $t=0$:
\begin{xmdsCode}
<filter>
  <vectors>main</vectors>
  <fourier_space>no</fourier_space>
  <![CDATA[
    phi *= exp(-(t*t));
  ]]>
</filter> 
\end{xmdsCode}

As in the \xmdsTag{vector} element the space in which the filter is to
be applied is specified in a \xmdsTag{fourier\_space} assignment.  \xmdsTag{moment\_group} and \xmdsTag{functions} elements may be used in the \xmdsTag{filter} element just as they are used in the \xmdsTag{integrate} element.

\section{The output element} 
\label{sec:outputElement}

\begin{xmdsCode}
<output>
  <filename>nlse.xsil</filename>
  <group>
    <sampling>
      <fourier_space> no </fourier_space>
      <lattice>       50 </lattice>
      <moments>pow_dens</moments>
      <![CDATA[
        pow_dens=~phi*phi;
      ]]>
    </sampling>
  </group>
</output>
\end{xmdsCode}

The \xmdsTag{output} element is compulsory and defines what will be
used as the output.

The \xmdsTag{filename} assignment is optional.  From version 1.2 \xmds
can produce binary output, as well as the default ascii format.  Both
output formats will produce an ascii file in XSIL
format~\cite{web:caltech_xsil}.  The default name for this file is the
\xmdsTag{simulation} \xmdsTag{name} appended with ``\ttt{.xsil}''.
With ascii output all of the data is contained within the XSIL file.
However, with binary output, the data is pointed to by the (ascii)
XSIL file, where the actual data is distributed among different files
labelled according to simulation name and output moment group.  In
general, the filenames are combined like this: \tsf{simulation name} +
\ttt{mg} + \tsf{moment group number} + \ttt{.dat}.  For instance, if
the simulation name is \ttt{nlse}, then the output binary file for the
first moment group will be called \ttt{nlsemg1.dat}.

In order that old output data remains meaningful, \xmds copies the
input simulation script to the name of the output file (overwriting it
if it already exists), and then the output data is inserted as
\xmdsTag{XSIL} elements before the closing \xmdsTag{/simulation}
tag. This way, the input parameters that generated the output data
remain known.

Next \xmds looks for the \xmdsTag{group} elements contained
within, of which there must be one or more. The \xmdsTag{group} element
contains moments that are sampled on a common output lattice in a
common Fourier space. If necessary the group may then be
post-processed after all segments have been performed, for example to
calculate moments with the propagation dimension itself in Fourier
space.

Thus \xmds looks for a compulsory \xmdsTag{sampling} element
within the \xmdsTag{group} element that defines how to sample for the
moments. It then looks for an optional \xmdsTag{post\_propagation}
element (within the same \xmdsTag{group} element) that defines any
post-processing required.

Consider firstly the \xmdsTag{sampling} element: Firstly, a
\xmdsTag{fourier\_space} assignment is compulsory, which contains as
many ``\ttt{yes}'' or ``\ttt{no}'' entries as there are transverse
dimensions, so that the field may be transformed accordingly prior to
sampling. In this example the field's ``\ttt{t}'' dimension remains
in normal space.

Next \xmds looks for a \xmdsTag{lattice} assignment containing the
lattice on which to sample on, for to sample the entire field may
produce much more data than is necessary to obtain a good plot. In
this assignment there should be as many integers as there are
transverse dimensions. Each integer should be one of the following
values (or else \xmds will exit with an appropriate error
message):

\begin{itemize}
\item{0:} A lattice integer of zero requests \xmds to integrate
the moments over this dimension. This will cause the output field to
no longer be a function of this transverse dimension.
\item{1:} A lattice integer of 1 requests \xmds to sample the
moments on a cross-sectional slice of this dimension, also causing the
output field to loose this transverse dimension. If this dimension is
in normal space then \xmds will extract the slice at the middle
lattice point (point number $N/2 + 1$ using integer division),
otherwise \xmds will extract the slice at the zero momentum
point, $k=0$.
\item A factor of the main field's corresponding lattice value which
is greater than one. In this case, if it is sampled in normal space
\xmds uses a coarser lattice that still spans the entire domain,
whereas if it is in Fourier space it samples a narrower window of the
k-space domain centered on the zero momentum point, $k=0$. You may of
course specify the same number of lattice points here as there are in
this dimension, in which case the output moments are simply mapped
point for point.
\end{itemize}

If there are no transverse dimensions, \xmds will not look for
either of the above two assignments. Next comes the \xmdsTag{moments}
assignment which lists the names of the moments. Any number of moments
are possible, provided there is at least one. Finally the code used
for calculating these moments is inserted as the content of the
\xmdsTag{sampling} element. The moments are complex variables and so
the code may be written accordingly.

Next, a \xmdsTag{post\_propagation} element may optionally follow the
\xmdsTag{sampling} element. In this element there must be a
\xmdsTag{fourier\_space} element with a ``\ttt{yes}'' or ``\ttt{no}''
entry for the propagation dimension, and as many more as there are
{\em remaining} transverse dimensions in the sampled output. In the
case of partially collapsed output the correspondence is
sequential. For example if the field had three transverse dimensions
``\ttt{x y z}'' and the``\ttt{y}'' dimension was collapsed by using
a \xmdsTag{lattice} assignment (within the \xmdsTag{sampling} element)
like ``\ttt{50 1 50}'', then a \xmdsTag{fourier\_space} assignment of
``\ttt{yes no yes}'' would mean transform to Fourier space in the
propagation dimension, normal space in the ``\ttt{x}'' dimension, and
Fourier space in the ``\ttt{z}'' dimension when evaluating the
post-processing moments. Finally there is another \xmdsTag{moments}
assignment (with new names), and more code in the content for
processing the moments already derived at the sampling level. As with
most other code blocks functions of dimensions may also be included,
but only of the dimensions that haven't been collapsed (plus the
propagation dimension).  Finally, it is only the real part of the
\xmdsTag{post\_propagation} moments that are written as output, or added
to the stochastic sums before they are processed.

In the absence of a \xmdsTag{post\_propagation} element the real parts
of the \xmdsTag{sampling} moments are used instead.