1<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01//EN">
2
3<!--Converted with LaTeX2HTML 2018.3 (Released July 19, 2018) -->
4<HTML lang="EN">
5<HEAD>
6<TITLE>synfast</TITLE>
7<META NAME="description" CONTENT="synfast">
8<META NAME="keywords" CONTENT="facilities">
9<META NAME="resource-type" CONTENT="document">
10<META NAME="distribution" CONTENT="global">
11
12<META HTTP-EQUIV="Content-Type" CONTENT="text/html; charset=utf-8">
13<META NAME="viewport" CONTENT="width=device-width, initial-scale=1.0">
14<META NAME="Generator" CONTENT="LaTeX2HTML v2018.3">
15   <link rel='apple-touch-icon' sizes='180x180' href='images/favicons/apple-touch-icon.png?v=2017'>
16   <link rel='icon' type='image/png' sizes='32x32' href='images/favicons/favicon-32x32.png?v=2017'>
17   <link rel='icon' type='image/png' sizes='16x16' href='images/favicons/favicon-16x16.png?v=2017'>
18   <link rel='manifest' href='images/favicons/manifest.json?v=2017'>
19   <link rel='mask-icon' href='images/favicons/safari-pinned-tab.svg?v=2017' color='#5bbad5'>
20   <link rel='shortcut icon' href='images/favicons/favicon.ico?v=2017'>
21   <meta name='apple-mobile-web-app-title' content='HEALPix'>
22   <meta name='application-name' content='HEALPix'>
23   <meta name='msapplication-config' content='images/favicons/browserconfig.xml?v=2017'>
24   <meta name='theme-color' content='#ffffff'>
25
26<LINK REL="STYLESHEET" HREF="facilities.css">
27
28<LINK REL="next" HREF="fac_ud_grade.htm">
29<LINK REL="previous" HREF="fac_smoothing.htm">
30<LINK REL="next" HREF="fac_ud_grade.htm">
31</HEAD>
32
33<body text="#000000" bgcolor="#FFFFFA">
34
35<DIV CLASS="navigation"><!--Navigation Panel-->
36<A
37 HREF="fac_smoothing.htm">
38<IMG WIDTH="63" HEIGHT="24" ALT="previous" SRC="prev.png"></A>
39<A
40 HREF="fac_HEALPix_F90_facilities.htm">
41<IMG WIDTH="26" HEIGHT="24" ALT="up" SRC="up.png"></A>
42<A
43 HREF="fac_ud_grade.htm">
44<IMG WIDTH="37" HEIGHT="24" ALT="next" SRC="next.png"></A>
45<A ID="tex2html133"
46  HREF="fac_TABLE_CONTENTS.htm">
47<IMG WIDTH="65" HEIGHT="24" ALT="contents" SRC="contents.png"></A>
48<BR>
49<B> Previous:</B> <A
50 HREF="fac_smoothing.htm">smoothing</A>
51
52<B>Up:</B> <A
53 HREF="fac_HEALPix_F90_facilities.htm">HEALPix/F90 facilities</A>
54
55<B> Next:</B> <A
56 HREF="fac_ud_grade.htm">ud_grade</A>
57<B> Top:</B> <a href="main.htm">Main Page</a></DIV>
58<!--End of Navigation Panel-->
59
60<H1><A ID="SECTION1400"></A>
61<A ID="fac:synfast"></A>
62<BR>
63
64</H1>
65
66<P>
67<b><font size=+6><FONT COLOR="#CC0000">synfast</FONT></font></b><hr>
68<H3>This program can be used to create  <b>HEALPix</b> maps (temperature only
69or temperature and polarisation)  computed as realisations
70of random Gaussian
71fields on a sphere characterized by the user provided
72theoretical power spectra,
73or as constrained realisations of such fields characterised by the user
74provided <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
75 ALT="$a_{\ell m}$"></SPAN> coefficients and/or power spectra.
76Total operation count scales as
77 <!-- MATH
78 ${\cal {O}}(N_{\mathrm{pix}}^{1/2}\ell_{\mathrm{max}}^2)$
79 -->
80<SPAN CLASS="MATH"><IMG STYLE="height: 3.21ex; vertical-align: -0.91ex; " SRC="fac_img26.png"
81 ALT="${\cal {O}}(N_{\mathrm{pix}}^{1/2}\ell_{\mathrm{max}}^2 )$"></SPAN> where <!-- MATH
82 $N_{\mathrm{pix}}$
83 -->
84<SPAN CLASS="MATH"><I>N</I><SUB>pix</SUB></SPAN> is the total number of pixels and
85<!-- MATH
86 $\ell_{\mathrm{max}}$
87 -->
88<SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img25.png"
89 ALT="$\ell_{\mathrm{max}}$"></SPAN> is the limiting spherical harmonics order.
90The map resolution, Gaussian beam FWHM,
91and random seed for the simulation can be selected by the user.
92Spherical harmonics are either generated using the recurrence relations
93during the execution of spectral synthesis, or  precomputed and read in
94before the synthesis is executed. The latter is no longer recommended since
95it provides no acceleration since the introduction of optimized algorithms
96in <b>HEALPix</b> v2.20. </H3>
97Location in HEALPix directory tree: <a href="https://sourceforge.net/p/healpix/code/1005/tree/trunk/src/f90/synfast/synfast.f90"><b>src/f90/synfast/synfast.f90</b></a>&nbsp;
98
99<P>
100<hr><h1>FORMAT </h1><blockquote><h3>%
101synfast [options] [parameter_file]
102</h3></blockquote>
103
104<P>
105<hr><H1>COMMAND LINE OPTIONS</H1>
106  <DL COMPACT><DT>
107<B><TT>-d</TT></B>
108<DD><DT>
109<B><TT>--double</TT></B>
110<DD>double precision mode (see
111<A HREF="fac_Using_HEALPix_Fortran_90_fa.htm#fac:subsec:ioprec">Notes on double/single precision modes</A>
112)
113    <DT>
114<B><TT>-s</TT></B>
115<DD><DT>
116<B><TT>--single</TT></B>
117<DD>single precision mode (default)
118  </DL>
119
120<P>
121<hr>
122<H1>QUALIFIERS</H1>
123
124  <DL COMPACT><DT>
125<B>infile = </B>
126<DD><A ID="fac:synfast:infile"></A>Defines the input power spectrum file,
127	(default= cl.fits). Note that <TT>infile</TT> is now optional :
128    <FONT COLOR="#CC0000">synfast</FONT> can run even if only <TT>almsfile</TT> is provided.
129    <DT>
130<B>outfile = </B>
131<DD><A ID="fac:synfast:outfile"></A>Defines the output (RING ordered) map file,
132(default= map.fits). Note that <TT>outfile</TT> is now optional: if it set to
133      `' (empty string),  mo map is synthesized but the <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
134 ALT="$a_{\ell m}$"></SPAN> generated can be output.
135    <DT>
136<B>outfile_alms = </B>
137<DD><A ID="fac:synfast:outfile_alms"></A>Defines the FITS file in which to output <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
138 ALT="$a_{\ell m}$"></SPAN> used
139      for the simulation (default= `')
140     <DT>
141<B>simul_type = </B>
142<DD><A ID="fac:synfast:simul_type"></A>Defines the simulation type, 1=temperature only (1 field),
143       2=temperature+polarisation (3 fields), 3=temperature and its first
144spatial derivatives (3 fields),
145       4=temperature and its first and second spatial derivatives (6 fields), 5=temperature
146       and polarisation, and their first derivatives (9 fields), 6=same as 5
147       plus the second derivatives of (T,Q,U) (18 fields).
148(default= 1).
149    <DT>
150<B>nsmax = </B>
151<DD><A ID="fac:synfast:nsmax"></A>Defines the resolution of the map.
152(default= 32)
153     <DT>
154<B>nlmax = </B>
155<DD><A ID="fac:synfast:nlmax"></A>Defines the maximum <SPAN CLASS="MATH"><IMG STYLE="height: 1.69ex; vertical-align: -0.10ex; " SRC="fac_img5.png"
156 ALT="$\ell$"></SPAN> value
157to be used in the simulation. WARNING: <!-- MATH
158 $\ell_{\mathrm{max}}$
159 -->
160<SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img25.png"
161 ALT="$\ell_{\mathrm{max}}$"></SPAN> can not exceed
162the value <SPAN CLASS="MATH"><IMG STYLE="height: 1.52ex; vertical-align: -0.10ex; " SRC="fac_img67.png"
163 ALT="$4\cdot$"></SPAN> <TT>nsmax</TT>, because the coefficients of the  average Fourier
164pixel window functions
165are precomputed and provided up to this limit.
166(default= 64)
167      <DT>
168<B>iseed = </B>
169<DD><A ID="fac:synfast:iseed"></A>Defines the random seed to be used
170for the generation of <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
171 ALT="$a_{\ell m}$"></SPAN>s from the power spectrum.
172(default= -1)
173    <DT>
174<B>fwhm_arcmin = </B>
175<DD><A ID="fac:synfast:fwhm_arcmin"></A>Defines the FWHM size in arcminutes
176of the simulated Gaussian beam.
177(default= 420.0)
178<DT>
179<B>beam_file = </B>
180<DD><A ID="fac:synfast:beam_file"></A>Defines the FITS file (see <A HREF="fac_Using_HEALPix_Fortran_90_fa.htm#fac:subsec:beamfiles">&rdquo;Beam window function files&rdquo; in introduction</A>) describing the
181    Legendre window
182    function of the circular beam to be used for the
183    simulation. If set to an existing file name, it will override the
184    <TT>fhwm_arcmin</TT> given above. (default=`')
185<DT>
186<B>almsfile = </B>
187<DD><A ID="fac:synfast:almsfile"></A>Defines the input filename for a file
188    containing <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
189 ALT="$a_{\ell m}$"></SPAN>s for constrained realisations.
190(default= `'). If <TT>apply_windows</TT> is <EM>false</EM>
191those <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
192 ALT="$a_{\ell m}$"></SPAN>s are used as they are, without being multiplied
193by the beam or pixel window function (with the assumption that they already have the
194    correct window functions). If <TT>apply_windows</TT> is <EM>true</EM>, the beam and
195    pixel window functions chosen above are applied to the constraining <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
196 ALT="$a_{\ell m}$"></SPAN> (with the
197    assumption that those are free of beam and pixel window function). The code
198    does not check the validity of these asumptions; if none is true, use the
199    <A HREF="fac_alteralm.htm#fac:alteralm">alteralm</A> facility to modify or remove
200    the window functions contained in the constraining <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
201 ALT="$a_{\ell m}$"></SPAN>.
202<DT>
203<B>apply_windows = </B>
204<DD><A ID="fac:synfast:apply_windows"></A>Determines how the constraining <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
205 ALT="$a_{\ell m}$"></SPAN> read from
206     <TT>almsfile</TT> are
207     treated with respect to window functions; see above for details.
208     y, yes, t, true, .true. and 1 are considered as <EM>true</EM>, while n, no, f,
209     false, .false. and 0 are considered as <EM>false</EM>, (default = .false.).
210<DT>
211<B>plmfile = </B>
212<DD><A ID="fac:synfast:plmfile"></A>Defines the input  filename for a file
213    containing  precomputed Legendre polynomials <SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img51.png"
214 ALT="$P_{\ell m}$"></SPAN>.
215(default= `')
216<DT>
217<B>windowfile = </B>
218<DD><A ID="fac:synfast:windowfile"></A>Defines the input filename  for the pixel
219    smoothing windows
220(default= pixel_window_n????.fits, see below).
221<DT>
222<B>winfiledir = </B>
223<DD><A ID="fac:synfast:winfiledir"></A>Defines the directory in which windowfile
224    is located (default: see
225<A HREF="fac_Using_HEALPix_Fortran_90_fa.htm#fac:subsec:defdir">Notes on default files and directories</A>).
226  </DL>
227
228<P>
229<hr>
230<H1>DESCRIPTION</H1>
231<blockquote>
232Synfast reads the power spectrum from a file in  ascii FITS
233format. This can contain either just the temperature power spectrum <SPAN CLASS="MATH"><IMG STYLE="height: 2.51ex; vertical-align: -0.67ex; " SRC="fac_img84.png"
234 ALT="$C^T_{\ell}$"></SPAN>s or
235temperature and polarisation power spectra: <SPAN CLASS="MATH"><IMG STYLE="height: 2.51ex; vertical-align: -0.67ex; " SRC="fac_img84.png"
236 ALT="$C^T_{\ell}$"></SPAN>, <SPAN CLASS="MATH"><IMG STYLE="height: 2.51ex; vertical-align: -0.67ex; " SRC="fac_img85.png"
237 ALT="$C^E_{\ell}$"></SPAN>, <SPAN CLASS="MATH"><IMG STYLE="height: 2.51ex; vertical-align: -0.67ex; " SRC="fac_img86.png"
238 ALT="$C^B_{\ell}$"></SPAN>
239and <!-- MATH
240 $C^{T\times E}_{\ell}$
241 -->
242<SPAN CLASS="MATH"><IMG STYLE="height: 2.62ex; vertical-align: -0.68ex; " SRC="fac_img87.png"
243 ALT="$C^{T\times E}_{\ell}$"></SPAN> (see <A HREF="#fac:synfast:note1">Note 1, below</A>). If <TT>simul_type = 2</TT> <FONT COLOR="#CC0000">synfast</FONT> generates
244Q and U maps as well as the temperature map. The output map(s)
245is (are) saved in a FITS file.
246The <SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img88.png"
247 ALT="$C_{\ell}$"></SPAN>s are used up to the specified
248<!-- MATH
249 $\ell_{\mathrm{max}}$
250 -->
251<SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img25.png"
252 ALT="$\ell_{\mathrm{max}}$"></SPAN>, which can not exceed <SPAN CLASS="MATH">4 x </SPAN> nsmax. If <TT>simul_type = 3</TT> or
253<TT>4</TT> the first derivatives of the temperature field or the first and second derivatives respectively
254are output as well as the temperature itself: <SPAN CLASS="MATH"><I>T</I>(<I>p</I>)</SPAN>, <!-- MATH
255 $\left({\partial T}/{\partial \theta}, {\partial T}/{\partial \phi}/\sin\theta \right)$
256 -->
257<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img89.png"
258 ALT="$\left({\partial T}/{\partial \theta}, {\partial T}/{\partial \phi}/\sin\theta \right)
259$"></SPAN>, <!-- MATH
260 $\left({\partial^2 T}/{\partial \theta^2}, {\partial^2 T}/{\partial
261  \theta\partial\phi}/\sin\theta,\right.$
262 -->
263<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img90.png"
264 ALT="$\left({\partial^2 T}/{\partial \theta^2}, {\partial^2 T}/{\partial
265\theta\partial\phi}/\sin\theta,\right. $"></SPAN>
266<!-- MATH
267 $\left.{\partial^2 T}/{\partial \phi^2}/\sin^2\theta \right)$
268 -->
269<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img91.png"
270 ALT="$\left.{\partial^2 T}/{\partial \phi^2}/\sin^2\theta \right) $"></SPAN>.
271If <TT>simul_type = 5</TT> or
272<TT>6</TT> the first derivatives of the (T,Q,U) fields or the first and second derivatives respectively
273are output as well as the field themself: <SPAN CLASS="MATH"><I>T</I>(<I>p</I>)</SPAN>,  <SPAN CLASS="MATH"><I>Q</I>(<I>p</I>)</SPAN>,  <SPAN CLASS="MATH"><I>U</I>(<I>p</I>)</SPAN>,
274<!-- MATH
275 $\left({\partial T}/{\partial \theta}, {\partial Q}/{\partial \theta}, {\partial
276  U}/{\partial \theta}; {\partial T}/{\partial \phi}/\sin\theta, \ldots \right)$
277 -->
278<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img92.png"
279 ALT="$\left({\partial T}/{\partial \theta}, {\partial Q}/{\partial \theta}, {\partial
280U}/{\partial \theta}; {\partial T}/{\partial \phi}/\sin\theta, \ldots \right)
281$"></SPAN>, <!-- MATH
282 $\left({\partial^2 T}/{\partial \theta^2},\ldots; {\partial^2 T}/{\partial
283  \theta\partial\phi}/\sin\theta,\ldots ;\right.$
284 -->
285<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img93.png"
286 ALT="$\left({\partial^2 T}/{\partial \theta^2},\ldots; {\partial^2 T}/{\partial
287\theta\partial\phi}/\sin\theta,\ldots ;\right. $"></SPAN>
288<!-- MATH
289 $\left.{\partial^2 T}/{\partial \phi^2}/\sin^2\theta \ldots \right)$
290 -->
291<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img94.png"
292 ALT="$\left.{\partial^2 T}/{\partial \phi^2}/\sin^2\theta \ldots \right) $"></SPAN>
293<br><br>The random sequence seed for generation of <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
294 ALT="$a_{\ell m}$"></SPAN> from the
295power spectrum should be non-zero integer. If 0 is provided, a seed is generated
296randomly by the code, based on the current date and time.
297The map can be convolved with a gaussian beam for which a beamsize can
298be specified, or for an arbitrary <EM>circular</EM> beam for which the
299Legendre transform is provided. The map is automatically convolved with a pixel window
300function. These are stored in FITS files  in
301the <TT>healpix/data</TT> directory. If <FONT COLOR="#CC0000">synfast</FONT> is not run in a directory
302which has these files, or from a directory which can reach these files
303by a `<TT>../data/</TT>' or `<TT>./data/</TT>' specification, the system
304variable <TT>HEALPIX</TT> is used to locate the main <b>HEALPix</b> directory
305and its <TT>data</TT> subdirectory is scanned. Failing this, the location of these
306files must be specified (using winfiledir). In the interactive mode this is
307requested only when necessary (see
308<A HREF="fac_Using_HEALPix_Fortran_90_fa.htm#fac:subsec:defdir">Notes on default directories</A>
309).
310<br><br>If some of the <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
311 ALT="$a_{\ell m}$"></SPAN> in the simulations are constrained eg. from observations, a FITS file
312with these <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
313 ALT="$a_{\ell m}$"></SPAN> can be read. This FITS file contains
314the <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
315 ALT="$a_{\ell m}$"></SPAN> for certain <SPAN CLASS="MATH"><IMG STYLE="height: 1.69ex; vertical-align: -0.10ex; " SRC="fac_img5.png"
316 ALT="$\ell$"></SPAN> and <SPAN CLASS="MATH"><I>m</I></SPAN> values
317and also the standard deviation for these <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
318 ALT="$a_{\ell m}$"></SPAN>. The sky
319realisation which <FONT COLOR="#CC0000">synfast</FONT> produces will be statistically consistent
320with the constraining <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
321 ALT="$a_{\ell m}$"></SPAN>.
322<br><br>The code can also be used
323to generate a set of <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
324 ALT="$a_{\ell m}$"></SPAN> matching the input power spectra, beam size and
325pixel size with or without actually synthesizing the map. Those <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
326 ALT="$a_{\ell m}$"></SPAN> can be
327used as an input (constraining <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
328 ALT="$a_{\ell m}$"></SPAN>) to another <FONT COLOR="#CC0000">synfast</FONT> run.
329<br><br>
330</blockquote>
331<blockquote>
332Spherical harmonics values in the synthesis are obtained from a
333recurrence on associated Legendre polynomials <!-- MATH
334 $P_{\ell m}(\theta)$
335 -->
336<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img52.png"
337 ALT="$P_{\ell m}(\theta)$"></SPAN>.
338This recurrence consumed most of the CPU time used by <FONT COLOR="#CC0000">synfast</FONT> up to version
3392.15. We have therefore included an option to load precomputed values for the
340<!-- MATH
341 $P_{\ell m}(\theta)$
342 -->
343<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img52.png"
344 ALT="$P_{\ell m}(\theta)$"></SPAN> from a file generated by the <b>HEALPix</b> facility
345<A HREF="fac_plmgen.htm#fac:plmgen">plmgen</A>. Since the introduction of accelerated spherical
346harmonic transforms in <b>HEALPix</b> v2.20, this feature is obsolete and should no
347longer be used.
348
349<P>
350<FONT COLOR="#CC0000">synfast</FONT> will issue a warning if the input FITS file for the power spectrum does
351not contain the keyword <TT>POLNORM</TT>. This keyword indicates that the convention
352used for polarization is consistent with CMBFAST (and consistent with <b>HEALPix</b>
3531.2). See the
354<A ID="tex2html25"
355  HREF="intro.htm"><b>HEALPix</b> Primer</A>
356for details on the
357polarization convention and the interface with CMBFAST. If the
358keyword is not found, <EM>no attempt will be made</EM> to renormalize the power
359spectrum.
360If the keyword is present, it will be inherited by the simulated map.
361
362<P>
363<A ID="fac:synfast:note1"></A><B>Note 1:</B> to allow the generation of maps (and <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
364 ALT="$a_{\ell m}$"></SPAN>) with <!-- MATH
365 $C^{T\times B}_{\ell} \ne 0$
366 -->
367<SPAN CLASS="MATH"><IMG STYLE="height: 2.62ex; vertical-align: -0.68ex; " SRC="fac_img95.png"
368 ALT="$C^{T\times B}_{\ell} \ne 0$"></SPAN> and/or <!-- MATH
369 $C^{E\times B}_{\ell} \ne 0$
370 -->
371<SPAN CLASS="MATH"><IMG STYLE="height: 2.62ex; vertical-align: -0.68ex; " SRC="fac_img96.png"
372 ALT="$C^{E\times B}_{\ell} \ne 0$"></SPAN>,
373see the subroutine <A HREF="./sub_create_alm.htm#sub:create_alm"><TT>create_alm</TT></A>.
374
375</blockquote>
376
377<P>
378<hr><H1>DATASETS</H1>
379<h3>The following datasets are involved in the <b><FONT COLOR="#CC0000">synfast</FONT></b>
380 processing.</h3>
381
382<TABLE CELLPADDING=3 BORDER="1">
383<TR><TH ALIGN="LEFT" VALIGN="TOP" WIDTH=150><SPAN  CLASS="textbf">Dataset</SPAN></TH>
384<TH ALIGN="LEFT" VALIGN="TOP" WIDTH=175><SPAN  CLASS="textbf">Description</SPAN></TH>
385</TR>
386<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=150>&nbsp;</TD>
387<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
388</TR>
389<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=150>
390/data/pixel_window_nxxxx.fits</TD>
391<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>Files containing pixel windows for
392                   various nsmax.</TD>
393</TR>
394<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=150>&nbsp;</TD>
395<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
396</TR>
397<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=150></TD>
398<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
399</TR>
400</TABLE>
401
402<P>
403<hr><H1>SUPPORT    </H1><H3>This section lists those routines and facilities (including those <i>external</i> to the <b>HEALPix</b> distribution) which can assist in the utilisation of <b><FONT COLOR="#CC0000">synfast</FONT></b>.</H3>
404<DL COMPACT><DT>
405<B><A HREF="./sub_generate_beam.htm#sub:generate_beam">generate_beam</A></B>
406<DD>This <b>HEALPix</b> Fortran
407subroutine generates or reads the <SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img21.png"
408 ALT="$B(\ell)$"></SPAN> window function used in <FONT COLOR="#CC0000">synfast</FONT>
409  <DT>
410<B><A HREF="fac_map2gif.htm#fac:map2gif">map2gif</A></B>
411<DD>This <b>HEALPix</b> Fortran facility can be used to visualise the
412  output map of <FONT COLOR="#CC0000">synfast</FONT>.
413  <DT>
414<B><A HREF="./idl_mollview.htm#idl:mollview">mollview</A></B>
415<DD>This <b>HEALPix</b> IDL facility can be used to visualise the
416  output map of <FONT COLOR="#CC0000">synfast</FONT>.
417  <DT>
418<B><A HREF="fac_alteralm.htm#fac:alteralm">alteralm</A></B>
419<DD>This <b>HEALPix</b> Fortran facility can be
420  used to implement the beam and pixel window functions on the constraining
421  <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
422 ALT="$a_{\ell m}$"></SPAN>s (<TT>almsfile</TT> file).
423  <DT>
424<B><A HREF="fac_anafast.htm#fac:anafast">anafast</A></B>
425<DD>This <b>HEALPix</b> Fortran facility can analyse a <b>HEALPix</b> map and
426     	       save the <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
427 ALT="$a_{\ell m}$"></SPAN> and <SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img88.png"
428 ALT="$C_{\ell}$"></SPAN>s to be read by <FONT COLOR="#CC0000">synfast</FONT>.
429  <DT>
430<B><A HREF="fac_plmgen.htm#fac:plmgen">plmgen</A></B>
431<DD>This <b>HEALPix</b> Fortran facility can be used to generate precomputed Legendre polynomials.
432
433<P>
434</DL>
435
436<P>
437<hr><H1>EXAMPLE # 1:</H1>
438<tt><TABLE CELLPADDING=3>
439<TR><TD ALIGN="LEFT">synfast</TD>
440<TD ALIGN="LEFT">&nbsp;</TD>
441</TR>
442</TABLE></tt>
443<blockquote>
444Synfast runs in interactive mode, self-explanatory.
445</blockquote>
446
447<P>
448
449<P>
450<hr><H1>EXAMPLE # 2:</H1>
451<tt><TABLE CELLPADDING=3>
452<TR><TD ALIGN="LEFT">synfast  filename</TD>
453<TD ALIGN="LEFT">&nbsp;</TD>
454</TR>
455</TABLE></tt>
456<blockquote>When 'filename' is present, <FONT COLOR="#CC0000">synfast</FONT> enters the non-interactive mode and parses
457its inputs from the file 'filename'. This has the following
458structure: the first entry is a qualifier which announces to the parser
459which input immediately follows. If this input is omitted in the
460input file, the parser assumes the default value.
461If the equality sign is omitted, then the parser ignores the entry.
462In this way comments may also be included in the file.
463In this example, the file contains the following qualifiers:
464<BR>
465<tt><A HREF="#fac:synfast:simul_type">simul_type</A>= 1</tt><br>
466<tt><A HREF="#fac:synfast:nsmax">nsmax</A>= 32</tt><br>
467<tt><A HREF="#fac:synfast:nlmax">nlmax</A>= 64</tt><br>
468<tt><A HREF="#fac:synfast:iseed">iseed</A>= -1</tt><br>
469<tt><A HREF="#fac:synfast:fwhm_arcmin">fwhm_arcmin</A>= 420.0</tt><br>
470<tt><A HREF="#fac:synfast:infile">infile</A>= cl.fits</tt><br>
471<tt><A HREF="#fac:synfast:outfile">outfile</A>= map.fits</tt><br>
472
473<P>
474Synfast reads in the <SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img88.png"
475 ALT="$C_{\ell}$"></SPAN> power spectrum in 'cl.fits' up to <SPAN CLASS="MATH"><IMG STYLE="height: 1.69ex; vertical-align: -0.10ex; " SRC="fac_img97.png"
476 ALT="$\ell=64$"></SPAN>, and
477produces the (RING ordered) map
478'map.fits' which has <!-- MATH
479 $N_{\mathrm{side}}=32$
480 -->
481<SPAN CLASS="MATH"><I>N</I><SUB>side</SUB>=32</SPAN>.
482The map is convolved with a beam of FWHM 420.0 arcminutes. The <!-- MATH
483 $\hyperref{iseed}{}{}{fac:synfast:iseed}=-1$
484 -->
485<SPAN CLASS="MATH"><IMG STYLE="height: 1.57ex; vertical-align: -0.29ex; " SRC="fac_img98.png"
486 ALT="$\hyperref{iseed}{}{}{fac:synfast:iseed}=-1$"></SPAN> sets
487the random seed for the realisation. A different <!-- MATH
488 $\hyperref{iseed}{}{}{fac:synfast:iseed}$
489 -->
490<SPAN CLASS="MATH"><IMG STYLE="vertical-align: 20px"
491 SRC="fac_img99.png"
492 ALT="$\hyperref{iseed}{}{}{fac:synfast:iseed}$"></SPAN> would have given a different
493realisation from the same power spectrum.
494
495<P>
496Since <BR>
497<tt><A HREF="#fac:synfast:outfile_alms">outfile_alms</A></tt><br>
498<tt><A HREF="#fac:synfast:almsfile">almsfile</A></tt><br>
499<tt><A HREF="#fac:synfast:apply_windows">apply_windows</A></tt><br>
500<tt><A HREF="#fac:synfast:plmfile">plmfile</A></tt><br>
501<tt><A HREF="#fac:synfast:beam_file">beam_file</A></tt><br>
502<tt><A HREF="#fac:synfast:windowfile">windowfile</A></tt><br>
503were omitted, they take their default values (empty strings).
504This means that no file for constrained realisation or precomputed
505Legendre polynomials are read, the <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
506 ALT="$a_{\ell m}$"></SPAN> generated in the process are not
507output, and <FONT COLOR="#CC0000">synfast</FONT> attempts to find the pixel
508window files in the default directories (see
509<A HREF="fac_Using_HEALPix_Fortran_90_fa.htm#fac:subsec:defdir">Notes on default files and directories</A>).
510</blockquote>
511
512<P>
513<hr><H1>RELEASE NOTES</H1><blockquote>
514  <DL COMPACT><DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
515 Initial release (<b>HEALPix</b> 0.90)
516    <DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
517 Optional non-interactive operation. Proper FITS file
518    support. Improved reccurence algorithm for <!-- MATH
519 $P_{\ell m}(\theta)$
520 -->
521<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img52.png"
522 ALT="$P_{\ell m}(\theta)$"></SPAN> which can compute to higher <SPAN CLASS="MATH"><IMG STYLE="height: 1.69ex; vertical-align: -0.10ex; " SRC="fac_img5.png"
523 ALT="$\ell$"></SPAN> values. Improved pixel windows averaged over
524    actual HEALPix pixels. New functionality: constrained realisations, precomputed
525    <SPAN CLASS="MATH"><IMG STYLE="height: 2.04ex; vertical-align: -0.45ex; " SRC="fac_img51.png"
526 ALT="$P_{\ell m}$"></SPAN>. (<b>HEALPix</b> 1.00)
527    <DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
528 New functionality: constrained realisations and pixel
529    windows are now available for polarization as well. Arbitrary
530    circular beams can be used. New parser (<b>HEALPix</b> 1.20)
531    <DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
532 New functionnality: the generated <SPAN CLASS="MATH"><IMG STYLE="height: 1.46ex; vertical-align: -0.45ex; " SRC="fac_img14.png"
533 ALT="$a_{\ell m}$"></SPAN> can be output, and the map
534    synthesis itself can be skipped. First and second derivatives of the
535    temperature field can be produced on demand.
536    <DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
537 New functionnality: First and second derivatives of the
538    <SPAN CLASS="MATH"><I>Q</I></SPAN> and <SPAN CLASS="MATH"><I>U</I></SPAN> Stokes field can be produced on demand.
539    <DT><DD><IMG WIDTH="14" HEIGHT="14" SRC="blueball.png" ALT="*">
540 Bug correction: corrected numerical errors on derivatives
541<!-- MATH
542 $\partial X/\partial\theta$
543 -->
544<SPAN CLASS="MATH"><IMG STYLE="height: 2.33ex; vertical-align: -0.68ex; " SRC="fac_img15.png"
545 ALT="$\partial X/\partial\theta$"></SPAN>,
546<!-- MATH
547 $\partial^2 X/(\partial\theta\partial\phi\sin\theta)$
548 -->
549<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img16.png"
550 ALT="$\partial^2 X/(\partial\theta\partial\phi\sin\theta)$"></SPAN>,
551<!-- MATH
552 $\partial^2 X/\partial \theta^2$
553 -->
554<SPAN CLASS="MATH"><IMG STYLE="height: 2.45ex; vertical-align: -0.68ex; " SRC="fac_img17.png"
555 ALT="$\partial^2 X/\partial \theta^2$"></SPAN>,
556for <SPAN CLASS="MATH"><I>X</I>=<I>Q</I>,<I>U</I></SPAN>. See <A HREF="fac_Appendix.htm#fac:sec:bug_synder">this appendix</A> for details.
557  (<b>HEALPix</b> 2.14)
558  </DL>
559</blockquote>
560
561<P>
562<hr><H1>MESSAGES</H1><h3>This section describes error messages generated by <b><FONT COLOR="#CC0000">synfast</FONT></b>
563</h3>
564
565<TABLE CELLPADDING=3 BORDER="1">
566<TR><TH ALIGN="LEFT" VALIGN="TOP" WIDTH=125><SPAN  CLASS="textbf">Message</SPAN></TH>
567<TH ALIGN="LEFT" VALIGN="TOP" WIDTH=50><SPAN  CLASS="textbf">Severity</SPAN></TH>
568<TH ALIGN="LEFT" VALIGN="TOP" WIDTH=175><SPAN  CLASS="textbf">Text</SPAN></TH>
569</TR>
570<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
571<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
572<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
573</TR>
574<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
575can not allocate memory for array xxx</TD>
576<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>Fatal</TD>
577<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>You do not have
578                   sufficient system resources to run this
579                   facility at the map resolution you required.
580  Try a lower map resolution.</TD>
581</TR>
582<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
583<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
584<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
585</TR>
586<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
587
588<P>
589this is not a binary table</TD>
590<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
591<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>the fitsfile you have specified is not
592of the proper format</TD>
593</TR>
594<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
595<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
596<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
597</TR>
598<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
599there are undefined values in the table!</TD>
600<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
601<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>the fitsfile you have specified is not
602of the proper format</TD>
603</TR>
604<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
605<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
606<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
607</TR>
608<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
609the header in xxx is too long</TD>
610<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
611<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>the fitsfile you have specified is not
612of the proper format</TD>
613</TR>
614<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
615<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
616<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
617</TR>
618<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
619XXX-keyword not found</TD>
620<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
621<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>the fitsfile you have specified is not
622of the proper format</TD>
623</TR>
624<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
625<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
626<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
627</TR>
628<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
629found xxx in the file, expected:yyyy</TD>
630<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
631<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>the specified fitsfile does not
632contain the proper amount of data.</TD>
633</TR>
634<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>&nbsp;</TD>
635<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
636<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
637</TR>
638<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=125>
639
640<P></TD>
641<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=50>&nbsp;</TD>
642<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=175>&nbsp;</TD>
643</TR>
644</TABLE>
645
646<P>
647
648<P>
649
650
651<P>
652
653<P>
654
655<DIV CLASS="navigation"><HR>
656<!--Navigation Panel-->
657<A
658 HREF="fac_smoothing.htm">
659<IMG WIDTH="63" HEIGHT="24" ALT="previous" SRC="prev.png"></A>
660<A
661 HREF="fac_HEALPix_F90_facilities.htm">
662<IMG WIDTH="26" HEIGHT="24" ALT="up" SRC="up.png"></A>
663<A
664 HREF="fac_ud_grade.htm">
665<IMG WIDTH="37" HEIGHT="24" ALT="next" SRC="next.png"></A>
666<A ID="tex2html133"
667  HREF="fac_TABLE_CONTENTS.htm">
668<IMG WIDTH="65" HEIGHT="24" ALT="contents" SRC="contents.png"></A>
669<BR>
670<B> Previous:</B> <A
671 HREF="fac_smoothing.htm">smoothing</A>
672
673<B>Up:</B> <A
674 HREF="fac_HEALPix_F90_facilities.htm">HEALPix/F90 facilities</A>
675
676<B> Next:</B> <A
677 HREF="fac_ud_grade.htm">ud_grade</A>
678<B> Top:</B> <a href="main.htm">Main Page</a></DIV>
679<!--End of Navigation Panel-->
680<ADDRESS>
681Version 3.50, 2018-12-10
682</ADDRESS>
683</BODY>
684</HTML>
685