1 | <!DOCTYPE html> |
---|
2 | <html lang="en"> |
---|
3 | <head> |
---|
4 | <meta charset="utf-8"> |
---|
5 | <meta http-equiv="X-UA-Compatible" content="IE=edge"> |
---|
6 | <meta name="viewport" content="width=device-width, initial-scale=1"> |
---|
7 | |
---|
8 | <meta name="description" content="Calculation of vertical velocity for FLEXPART"> |
---|
9 | |
---|
10 | <meta name="author" content="Leopold Haimberger" > |
---|
11 | <link rel="icon" href="../favicon.png"> |
---|
12 | |
---|
13 | <title>phgrreal.f – Flex_extract: Calculation of etadot</title> |
---|
14 | |
---|
15 | <link href="../css/bootstrap.min.css" rel="stylesheet"> |
---|
16 | <link href="../css/pygments.css" rel="stylesheet"> |
---|
17 | <link href="../css/font-awesome.min.css" rel="stylesheet"> |
---|
18 | <link href="../css/local.css" rel="stylesheet"> |
---|
19 | |
---|
20 | |
---|
21 | <!-- HTML5 shim and Respond.js for IE8 support of HTML5 elements and media queries --> |
---|
22 | <!--[if lt IE 9]> |
---|
23 | <script src="https://oss.maxcdn.com/html5shiv/3.7.2/html5shiv.min.js"></script> |
---|
24 | <script src="https://oss.maxcdn.com/respond/1.4.2/respond.min.js"></script> |
---|
25 | <![endif]--> |
---|
26 | |
---|
27 | <script src="../js/jquery-2.1.3.min.js"></script> |
---|
28 | <script src="../js/svg-pan-zoom.min.js"></script> |
---|
29 | |
---|
30 | |
---|
31 | <script src="../tipuesearch/tipuesearch_content.js"></script> |
---|
32 | <link href="../tipuesearch/tipuesearch.css" rel="stylesheet"> |
---|
33 | <script src="../tipuesearch/tipuesearch_set.js"></script> |
---|
34 | <script src="../tipuesearch/tipuesearch.js"></script> |
---|
35 | |
---|
36 | |
---|
37 | </head> |
---|
38 | |
---|
39 | <body> |
---|
40 | |
---|
41 | <!-- Fixed navbar --> |
---|
42 | <nav class="navbar navbar-inverse navbar-fixed-top"> |
---|
43 | <div class="container"> |
---|
44 | <div class="navbar-header"> |
---|
45 | <button type="button" class="navbar-toggle collapsed" data-toggle="collapse" data-target="#navbar" aria-expanded="false" aria-controls="navbar"> |
---|
46 | <span class="sr-only">Toggle navigation</span> |
---|
47 | <span class="icon-bar"></span> |
---|
48 | <span class="icon-bar"></span> |
---|
49 | <span class="icon-bar"></span> |
---|
50 | </button> |
---|
51 | <a class="navbar-brand" href="../index.html">Flex_extract: Calculation of etadot <small>7.1</small></a> |
---|
52 | </div> |
---|
53 | <div id="navbar" class="navbar-collapse collapse"> |
---|
54 | <ul class="nav navbar-nav"> |
---|
55 | |
---|
56 | <li class="dropdown hidden-xs visible-sm visible-md hidden-lg"> |
---|
57 | <a href="#" class="dropdown-toggle" |
---|
58 | data-toggle="dropdown" role="button" |
---|
59 | aria-haspopup="true" |
---|
60 | aria-expanded="false">Contents <span class="caret"></span></a> |
---|
61 | <ul class="dropdown-menu"> |
---|
62 | |
---|
63 | <li><a href="../lists/files.html">Source Files</a></li> |
---|
64 | |
---|
65 | |
---|
66 | <li><a href="../lists/modules.html">Modules</a></li> |
---|
67 | |
---|
68 | |
---|
69 | |
---|
70 | <li><a href="../lists/procedures.html">Procedures</a></li> |
---|
71 | |
---|
72 | |
---|
73 | |
---|
74 | <li><a href="../program/preconvert.html">Program</a></li> |
---|
75 | |
---|
76 | </ul> |
---|
77 | </li> |
---|
78 | |
---|
79 | <li class="visible-xs hidden-sm visible-lg"><a href="../lists/files.html">Source Files</a></li> |
---|
80 | |
---|
81 | |
---|
82 | <li class="visible-xs hidden-sm visible-lg"><a href="../lists/modules.html">Modules</a></li> |
---|
83 | |
---|
84 | |
---|
85 | |
---|
86 | <li class="visible-xs hidden-sm visible-lg"><a href="../lists/procedures.html">Procedures</a></li> |
---|
87 | |
---|
88 | |
---|
89 | |
---|
90 | <li class="visible-xs hidden-sm visible-lg"><a href="../program/preconvert.html">Program</a></li> |
---|
91 | |
---|
92 | </ul> |
---|
93 | |
---|
94 | <form action="../search.html" class="navbar-form navbar-right" role="search"> |
---|
95 | <div class="form-group"> |
---|
96 | <input type="text" class="form-control" placeholder="Search" name="q" id="tipue_search_input" autocomplete="off" required> |
---|
97 | </div> |
---|
98 | <!-- |
---|
99 | <button type="submit" class="btn btn-default">Submit</button> |
---|
100 | --> |
---|
101 | </form> |
---|
102 | |
---|
103 | </div><!--/.nav-collapse --> |
---|
104 | </div> |
---|
105 | </nav> |
---|
106 | |
---|
107 | <div class="container"> |
---|
108 | |
---|
109 | |
---|
110 | <div class="row"> |
---|
111 | <h1>phgrreal.f |
---|
112 | <small>Source File</small> |
---|
113 | |
---|
114 | </h1> |
---|
115 | |
---|
116 | <div class="row"> |
---|
117 | <div class="col-lg-12"> |
---|
118 | <div class="well well-sm"> |
---|
119 | <ul class="list-inline" style="margin-bottom:0px;display:inline"> |
---|
120 | |
---|
121 | |
---|
122 | |
---|
123 | |
---|
124 | |
---|
125 | |
---|
126 | <li><i class="fa fa-list-ol"></i> |
---|
127 | <a data-toggle="tooltip" |
---|
128 | data-placement="bottom" data-html="true" |
---|
129 | title="29.9% of total for source files.">397 statements</a> |
---|
130 | </li> |
---|
131 | |
---|
132 | |
---|
133 | <li><i class="fa fa-code"></i><a href="../src/phgrreal.f"> Source File</a></li> |
---|
134 | |
---|
135 | </ul> |
---|
136 | <ol class="breadcrumb in-well text-right"> |
---|
137 | |
---|
138 | <li class="active">phgrreal.f</li> |
---|
139 | </ol> |
---|
140 | </div> |
---|
141 | </div> |
---|
142 | </div> |
---|
143 | <script> |
---|
144 | $(function () { |
---|
145 | $('[data-toggle="tooltip"]').tooltip() |
---|
146 | }) |
---|
147 | </script> |
---|
148 | |
---|
149 | </div> |
---|
150 | <div class="row"> |
---|
151 | <div class="col-md-3 hidden-xs hidden-sm visible-md visible-lg"> |
---|
152 | |
---|
153 | <div id="sidebar"> |
---|
154 | |
---|
155 | |
---|
156 | |
---|
157 | |
---|
158 | |
---|
159 | |
---|
160 | |
---|
161 | <div class="panel panel-primary"> |
---|
162 | <div class="panel-heading text-left"><h3 class="panel-title"><a data-toggle="collapse" href="#mods-0">Modules</a></h3></div> |
---|
163 | <div id="mods-0" class="panel-collapse collapse"> |
---|
164 | <div class="list-group"> |
---|
165 | |
---|
166 | <a class="list-group-item" href="../module/phtogr.html">PHTOGR</a> |
---|
167 | |
---|
168 | </div> |
---|
169 | </div> |
---|
170 | </div> |
---|
171 | |
---|
172 | |
---|
173 | |
---|
174 | |
---|
175 | |
---|
176 | |
---|
177 | |
---|
178 | |
---|
179 | |
---|
180 | |
---|
181 | |
---|
182 | |
---|
183 | |
---|
184 | |
---|
185 | |
---|
186 | <div class="panel panel-primary"> |
---|
187 | <div class="panel-heading text-left"><h3 class="panel-title">Source Code</h3></div> |
---|
188 | <div class="list-group"> |
---|
189 | <a class="list-group-item" href="../sourcefile/phgrreal.f.html#src">phgrreal.f</a> |
---|
190 | </div> |
---|
191 | </div> |
---|
192 | |
---|
193 | |
---|
194 | <hr> |
---|
195 | |
---|
196 | |
---|
197 | <div class="panel panel-default"> |
---|
198 | <div class="panel-heading text-left"><h3 class="panel-title"><a data-toggle="collapse" href="#allfiles-0">All Source Files</a></h3></div> |
---|
199 | <div id="allfiles-0" class="panel-collapse collapse"> |
---|
200 | <div class="list-group"> |
---|
201 | |
---|
202 | <a class="list-group-item" href="../sourcefile/grphreal.f.html">grphreal.f</a> |
---|
203 | |
---|
204 | <a class="list-group-item" href="../sourcefile/phgrreal.f.html">phgrreal.f</a> |
---|
205 | |
---|
206 | <a class="list-group-item" href="../sourcefile/preconvert.f90.html">preconvert.f90</a> |
---|
207 | |
---|
208 | <a class="list-group-item" href="../sourcefile/rwgrib2.f90.html">rwGRIB2.f90</a> |
---|
209 | |
---|
210 | </div> |
---|
211 | </div> |
---|
212 | </div> |
---|
213 | |
---|
214 | |
---|
215 | </div> |
---|
216 | |
---|
217 | </div> |
---|
218 | <div class="col-md-9" id='text'> |
---|
219 | |
---|
220 | |
---|
221 | |
---|
222 | <h3>Files Dependent On This One</h3> |
---|
223 | |
---|
224 | <div class="depgraph"><?xml version="1.0" encoding="UTF-8" standalone="no"?> |
---|
225 | <!DOCTYPE svg PUBLIC "-//W3C//DTD SVG 1.1//EN" |
---|
226 | "http://www.w3.org/Graphics/SVG/1.1/DTD/svg11.dtd"> |
---|
227 | <!-- Generated by graphviz version 2.38.0 (20140413.2041) |
---|
228 | --> |
---|
229 | <!-- Title: sourcefile~~phgrreal.f~~AfferentGraph Pages: 1 --> |
---|
230 | <svg id="sourcefilephgrrealfAfferentGraph" width="278pt" height="52pt" |
---|
231 | viewBox="0.00 0.00 278.00 52.00" xmlns="http://www.w3.org/2000/svg" xmlns:xlink="http://www.w3.org/1999/xlink"> |
---|
232 | <g id="sourcefile~~phgrreal.f~~AfferentGraph" class="graph" transform="scale(1 1) rotate(0) translate(4 48)"> |
---|
233 | <title>sourcefile~~phgrreal.f~~AfferentGraph</title> |
---|
234 | <polygon fill="white" stroke="none" points="-4,4 -4,-48 274,-48 274,4 -4,4"/> |
---|
235 | <!-- sourcefile~phgrreal.f --> |
---|
236 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_node1" class="node"><title>sourcefile~phgrreal.f</title> |
---|
237 | <polygon fill="none" stroke="black" points="58,-24 3.55271e-15,-24 3.55271e-15,-0 58,-0 58,-24"/> |
---|
238 | <text text-anchor="middle" x="29" y="-9.6" font-family="Helvetica,sans-Serif" font-size="10.50">phgrreal.f</text> |
---|
239 | </g> |
---|
240 | <!-- sourcefile~grphreal.f --> |
---|
241 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_node2" class="node"><title>sourcefile~grphreal.f</title> |
---|
242 | <g id="a_sourcefile~~phgrreal.f~~AfferentGraph_node2"><a xlink:href="../sourcefile/grphreal.f.html" xlink:title="grphreal.f"> |
---|
243 | <polygon fill="#f0ad4e" stroke="#f0ad4e" points="152,-44 94,-44 94,-20 152,-20 152,-44"/> |
---|
244 | <text text-anchor="middle" x="123" y="-29.6" font-family="Helvetica,sans-Serif" font-size="10.50" fill="white">grphreal.f</text> |
---|
245 | </a> |
---|
246 | </g> |
---|
247 | </g> |
---|
248 | <!-- sourcefile~phgrreal.f->sourcefile~grphreal.f --> |
---|
249 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_edge2" class="edge"><title>sourcefile~phgrreal.f->sourcefile~grphreal.f</title> |
---|
250 | <path fill="none" stroke="#000000" stroke-dasharray="5,2" d="M58.1029,-18.1093C66.2308,-19.8763 75.2419,-21.8352 83.872,-23.7113"/> |
---|
251 | <polygon fill="#000000" stroke="#000000" points="83.2368,-27.1549 93.7521,-25.8592 84.7239,-20.3147 83.2368,-27.1549"/> |
---|
252 | </g> |
---|
253 | <!-- sourcefile~preconvert.f90 --> |
---|
254 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_node3" class="node"><title>sourcefile~preconvert.f90</title> |
---|
255 | <g id="a_sourcefile~~phgrreal.f~~AfferentGraph_node3"><a xlink:href="../sourcefile/preconvert.f90.html" xlink:title="preconvert.f90"> |
---|
256 | <polygon fill="#f0ad4e" stroke="#f0ad4e" points="270,-24 188,-24 188,-0 270,-0 270,-24"/> |
---|
257 | <text text-anchor="middle" x="229" y="-9.6" font-family="Helvetica,sans-Serif" font-size="10.50" fill="white">preconvert.f90</text> |
---|
258 | </a> |
---|
259 | </g> |
---|
260 | </g> |
---|
261 | <!-- sourcefile~phgrreal.f->sourcefile~preconvert.f90 --> |
---|
262 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_edge3" class="edge"><title>sourcefile~phgrreal.f->sourcefile~preconvert.f90</title> |
---|
263 | <path fill="none" stroke="#000000" stroke-dasharray="5,2" d="M58.0369,-11.4671C69.224,-11.2793 82.2211,-11.092 94,-11 119.777,-10.7986 126.223,-10.8304 152,-11 160.255,-11.0543 169.015,-11.1418 177.539,-11.243"/> |
---|
264 | <polygon fill="#000000" stroke="#000000" points="177.726,-14.7456 187.77,-11.372 177.815,-7.74611 177.726,-14.7456"/> |
---|
265 | </g> |
---|
266 | <!-- sourcefile~grphreal.f->sourcefile~preconvert.f90 --> |
---|
267 | <g id="sourcefile~~phgrreal.f~~AfferentGraph_edge1" class="edge"><title>sourcefile~grphreal.f->sourcefile~preconvert.f90</title> |
---|
268 | <path fill="none" stroke="#000000" stroke-dasharray="5,2" d="M152.144,-26.5877C160.099,-25.058 169.004,-23.3453 177.839,-21.6463"/> |
---|
269 | <polygon fill="#000000" stroke="#000000" points="178.57,-25.07 187.729,-19.7444 177.248,-18.1959 178.57,-25.07"/> |
---|
270 | </g> |
---|
271 | </g> |
---|
272 | </svg> |
---|
273 | </div> |
---|
274 | <div><a type="button" class="graph-help" data-toggle="modal" href="#graph-help-text">Help</a></div> |
---|
275 | <div class="modal fade" id="graph-help-text" tabindex="-1" role="dialog"> |
---|
276 | <div class="modal-dialog modal-lg" role="document"> |
---|
277 | <div class="modal-content"> |
---|
278 | <div class="modal-header"> |
---|
279 | <button type="button" class="close" data-dismiss="modal" aria-label="Close"><span aria-hidden="true">×</span></button> |
---|
280 | <h4 class="modal-title" id="-graph-help-label">Graph Key</h4> |
---|
281 | </div> |
---|
282 | <div class="modal-body"> |
---|
283 | |
---|
284 | <p>Nodes of different colours represent the following: </p> |
---|
285 | <?xml version="1.0" encoding="UTF-8" standalone="no"?> |
---|
286 | <!DOCTYPE svg PUBLIC "-//W3C//DTD SVG 1.1//EN" |
---|
287 | "http://www.w3.org/Graphics/SVG/1.1/DTD/svg11.dtd"> |
---|
288 | <!-- Generated by graphviz version 2.38.0 (20140413.2041) |
---|
289 | --> |
---|
290 | <!-- Title: Graph Key Pages: 1 --> |
---|
291 | <svg width="190pt" height="32pt" |
---|
292 | viewBox="0.00 0.00 190.00 32.00" xmlns="http://www.w3.org/2000/svg" xmlns:xlink="http://www.w3.org/1999/xlink"> |
---|
293 | <g id="graph0" class="graph" transform="scale(1 1) rotate(0) translate(4 28)"> |
---|
294 | <title>Graph Key</title> |
---|
295 | <polygon fill="white" stroke="none" points="-4,4 -4,-28 186,-28 186,4 -4,4"/> |
---|
296 | <!-- Source File --> |
---|
297 | <g id="node1" class="node"><title>Source File</title> |
---|
298 | <polygon fill="#f0ad4e" stroke="#f0ad4e" points="67,-24 0,-24 0,-0 67,-0 67,-24"/> |
---|
299 | <text text-anchor="middle" x="33.5" y="-9.6" font-family="Helvetica,sans-Serif" font-size="10.50" fill="white">Source File</text> |
---|
300 | </g> |
---|
301 | <!-- This Page's Entity --> |
---|
302 | <g id="node2" class="node"><title>This Page's Entity</title> |
---|
303 | <polygon fill="none" stroke="black" points="182,-24 85,-24 85,-0 182,-0 182,-24"/> |
---|
304 | <text text-anchor="middle" x="133.5" y="-9.6" font-family="Helvetica,sans-Serif" font-size="10.50">This Page's Entity</text> |
---|
305 | </g> |
---|
306 | </g> |
---|
307 | </svg> |
---|
308 | |
---|
309 | |
---|
310 | <p>Solid arrows point from a file to a file which depends upon it. A file |
---|
311 | is dependent upon another if the latter must be compiled before the former |
---|
312 | can be. |
---|
313 | </p> |
---|
314 | |
---|
315 | </div> |
---|
316 | </div> |
---|
317 | </div> |
---|
318 | </div> |
---|
319 | |
---|
320 | |
---|
321 | <br> |
---|
322 | |
---|
323 | <section class="visible-xs visible-sm hidden-md"> |
---|
324 | |
---|
325 | |
---|
326 | |
---|
327 | |
---|
328 | |
---|
329 | |
---|
330 | |
---|
331 | <div class="panel panel-primary"> |
---|
332 | <div class="panel-heading text-left"><h3 class="panel-title"><a data-toggle="collapse" href="#mods-1">Modules</a></h3></div> |
---|
333 | <div id="mods-1" class="panel-collapse collapse"> |
---|
334 | <div class="list-group"> |
---|
335 | |
---|
336 | <a class="list-group-item" href="../module/phtogr.html">PHTOGR</a> |
---|
337 | |
---|
338 | </div> |
---|
339 | </div> |
---|
340 | </div> |
---|
341 | |
---|
342 | |
---|
343 | |
---|
344 | |
---|
345 | |
---|
346 | |
---|
347 | |
---|
348 | |
---|
349 | |
---|
350 | |
---|
351 | |
---|
352 | |
---|
353 | |
---|
354 | |
---|
355 | |
---|
356 | <div class="panel panel-primary"> |
---|
357 | <div class="panel-heading text-left"><h3 class="panel-title">Source Code</h3></div> |
---|
358 | <div class="list-group"> |
---|
359 | <a class="list-group-item" href="../sourcefile/phgrreal.f.html#src">phgrreal.f</a> |
---|
360 | </div> |
---|
361 | </div> |
---|
362 | |
---|
363 | |
---|
364 | </section> |
---|
365 | <br class="visible-xs visible-sm hidden-md"> |
---|
366 | |
---|
367 | <section> |
---|
368 | <h2><span class="anchor" id="src"></span>Source Code</h2> |
---|
369 | <div class="hl"><pre><span></span><a name="ln-1"></a><span class="nl"> </span> <span class="k">MODULE </span><span class="n">PHTOGR</span> |
---|
370 | <a name="ln-2"></a> |
---|
371 | <a name="ln-3"></a><span class="nl"> </span> <span class="kt">INTEGER</span><span class="p">,</span> <span class="k">PARAMETER</span> <span class="kd">::</span> <span class="n">MAXAUF</span><span class="o">=</span><span class="mi">36000</span> |
---|
372 | <a name="ln-4"></a> |
---|
373 | <a name="ln-5"></a><span class="nl"> </span> <span class="k">CONTAINS</span> |
---|
374 | <a name="ln-6"></a> |
---|
375 | <a name="ln-7"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PHGR213</span><span class="p">(</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">MLAT</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span> |
---|
376 | <a name="ln-8"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
377 | <a name="ln-9"></a> |
---|
378 | <a name="ln-10"></a><span class="c">C DIE ROUTINE F]HRT EINE TRANSFORMATION EINER</span> |
---|
379 | <a name="ln-11"></a><span class="c">C FELDVARIABLEN VOM PHASENRAUM IN DEN PHYSIKALISCHEN</span> |
---|
380 | <a name="ln-12"></a><span class="c">C RAUM AUF DAS REDUZIERTE GAUSS'SCHE GITTER DURCH</span> |
---|
381 | <a name="ln-13"></a><span class="c">C</span> |
---|
382 | <a name="ln-14"></a><span class="c">C CXMN = SPEKTRALKOEFFIZIENTEN IN DER REIHENFOLGE</span> |
---|
383 | <a name="ln-15"></a><span class="c">C CX00,CX01,CX11,CX02,....CXMNAUFMNAUF</span> |
---|
384 | <a name="ln-16"></a><span class="c">C FELD = FELD DER METEOROLOGISCHEN VARIABLEN</span> |
---|
385 | <a name="ln-17"></a><span class="c">C WSAVE = Working Array fuer Fouriertransformation</span> |
---|
386 | <a name="ln-18"></a><span class="c">C Z = LEGENDREFUNKTIONSWERTE</span> |
---|
387 | <a name="ln-19"></a><span class="c">C</span> |
---|
388 | <a name="ln-20"></a><span class="c">C MNAUF ANZAHL DER FOURIERKOEFFIZIENTEN</span> |
---|
389 | <a name="ln-21"></a><span class="c">C MAXL ANZAHL DER FUER DAS GITTER BENUTZTEN LAENGEN</span> |
---|
390 | <a name="ln-22"></a><span class="c">C MAXB ANZAHL DER FUER DAS GITTER BENOETIGTEN BREITEN</span> |
---|
391 | <a name="ln-23"></a><span class="c">C MLEVEL ANZAHL DER LEVELS, DIE TRANSFORMIERT WERDEN</span> |
---|
392 | <a name="ln-24"></a><span class="c">C</span> |
---|
393 | <a name="ln-25"></a><span class="nl"> </span> <span class="k">IMPLICIT NONE</span> |
---|
394 | <a name="ln-26"></a> |
---|
395 | <a name="ln-27"></a><span class="c">C Anzahl der Gitterpunkte auf jedem Breitenkreis</span> |
---|
396 | <a name="ln-28"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">MLAT</span><span class="p">(</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
397 | <a name="ln-29"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">K</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
398 | <a name="ln-30"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">IND</span><span class="p">(</span><span class="n">MAXB</span><span class="p">)</span> |
---|
399 | <a name="ln-31"></a><span class="nl"> </span> |
---|
400 | <a name="ln-32"></a> |
---|
401 | <a name="ln-33"></a><span class="c">C FELD DER LEGENDREPOLYNOME FUER EINE BREITE</span> |
---|
402 | <a name="ln-34"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
403 | <a name="ln-35"></a> |
---|
404 | <a name="ln-36"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
405 | <a name="ln-37"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
406 | <a name="ln-38"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">WSAVE</span><span class="p">(</span><span class="mi">8</span><span class="o">*</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">15</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
407 | <a name="ln-39"></a><span class="nl"> </span> <span class="kt">INTEGER</span> <span class="kd">::</span> <span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span><span class="n">MAXB</span><span class="p">)</span> |
---|
408 | <a name="ln-40"></a><span class="nl"> </span> |
---|
409 | <a name="ln-41"></a><span class="nl"> </span> <span class="n">IND</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="mi">0</span> |
---|
410 | <a name="ln-42"></a><span class="nl"> </span><span class="gs"> </span><span class="k">DO </span><span class="mi">7</span> <span class="n">K</span><span class="o">=</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span> |
---|
411 | <a name="ln-43"></a><span class="nl"> </span> <span class="n">IND</span><span class="p">(</span><span class="n">K</span><span class="p">)</span><span class="o">=</span><span class="n">IND</span><span class="p">(</span><span class="n">K</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span><span class="o">+</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
412 | <a name="ln-44"></a><span class="nl">7 C</span><span class="gs">O</span><span class="n">NTINUE</span> |
---|
413 | <a name="ln-45"></a> |
---|
414 | <a name="ln-46"></a><span class="c">!$OMP PARALLEL DO SCHEDULE(DYNAMIC)</span> |
---|
415 | <a name="ln-47"></a><span class="nl"> </span><span class="gs"> </span><span class="k">DO </span><span class="mi">17</span> <span class="n">K</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span> |
---|
416 | <a name="ln-48"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">PHSYM</span><span class="p">(</span><span class="n">K</span><span class="p">,</span><span class="n">IND</span><span class="p">,</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MLAT</span><span class="p">,</span> |
---|
417 | <a name="ln-49"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
418 | <a name="ln-50"></a> |
---|
419 | <a name="ln-51"></a><span class="nl">17 </span><span class="gs">C</span><span class="n">ONTINUE</span> |
---|
420 | <a name="ln-52"></a><span class="c">!$OMP END PARALLEL DO</span> |
---|
421 | <a name="ln-53"></a> |
---|
422 | <a name="ln-54"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
423 | <a name="ln-55"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PHGR213</span> |
---|
424 | <a name="ln-56"></a><span class="c">C</span> |
---|
425 | <a name="ln-57"></a><span class="c">C</span> |
---|
426 | <a name="ln-58"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PHSYM</span><span class="p">(</span><span class="n">K</span><span class="p">,</span><span class="n">IND</span><span class="p">,</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MLAT</span><span class="p">,</span> |
---|
427 | <a name="ln-59"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
428 | <a name="ln-60"></a> |
---|
429 | <a name="ln-61"></a><span class="nl"> </span> <span class="k">IMPLICIT NONE</span> |
---|
430 | <a name="ln-62"></a> |
---|
431 | <a name="ln-63"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">MLAT</span><span class="p">(</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
432 | <a name="ln-64"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">K</span><span class="p">,</span><span class="n">L</span><span class="p">,</span><span class="n">I</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">LLS</span><span class="p">,</span><span class="n">LLPS</span><span class="p">,</span><span class="n">LL</span><span class="p">,</span><span class="n">LLP</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
433 | <a name="ln-65"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">IND</span><span class="p">(</span><span class="n">MAXB</span><span class="p">)</span> |
---|
434 | <a name="ln-66"></a><span class="nl"> </span> <span class="kt">INTEGER</span> <span class="kd">::</span> <span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span><span class="n">MAXB</span><span class="p">)</span> |
---|
435 | <a name="ln-67"></a><span class="nl"> </span> |
---|
436 | <a name="ln-68"></a> |
---|
437 | <a name="ln-69"></a><span class="c">C FELD DER FOURIERKOEFFIZIENTEN</span> |
---|
438 | <a name="ln-70"></a><span class="nl"> </span> <span class="kt">REAL</span> <span class="kd">::</span> <span class="n">CXMS</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">),</span><span class="n">CXMA</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
439 | <a name="ln-71"></a> |
---|
440 | <a name="ln-72"></a><span class="c">C FELD DER LEGENDREPOLYNOME FUER EINE BREITE</span> |
---|
441 | <a name="ln-73"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
442 | <a name="ln-74"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">ACR</span><span class="p">,</span><span class="n">ACI</span><span class="p">,</span><span class="n">SCR</span><span class="p">,</span><span class="n">SCI</span> |
---|
443 | <a name="ln-75"></a> |
---|
444 | <a name="ln-76"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
445 | <a name="ln-77"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
446 | <a name="ln-78"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">WSAVE</span><span class="p">(</span><span class="mi">8</span><span class="o">*</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">15</span><span class="p">,</span><span class="n">MAXB</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
447 | <a name="ln-79"></a> |
---|
448 | <a name="ln-80"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">6</span> <span class="n">L</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span> |
---|
449 | <a name="ln-81"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="mi">0</span> |
---|
450 | <a name="ln-82"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="mi">0</span> |
---|
451 | <a name="ln-83"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">1</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
452 | <a name="ln-84"></a><span class="nl"> </span> <span class="n">SCR</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
453 | <a name="ln-85"></a><span class="nl"> </span> <span class="n">SCI</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
454 | <a name="ln-86"></a><span class="nl"> </span> <span class="n">ACR</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
455 | <a name="ln-87"></a><span class="nl"> </span> <span class="n">ACI</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
456 | <a name="ln-88"></a><span class="nl"> </span> <span class="n">LLS</span><span class="o">=</span><span class="n">LL</span> |
---|
457 | <a name="ln-89"></a><span class="nl"> </span> <span class="n">LLPS</span><span class="o">=</span><span class="n">LLP</span> |
---|
458 | <a name="ln-90"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mf">1.</span><span class="n">LT</span><span class="p">.</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">))</span> <span class="k">THEN</span> |
---|
459 | <a name="ln-91"></a><span class="c">C Innerste Schleife aufgespalten um if-Abfrage zu sparen</span> |
---|
460 | <a name="ln-92"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">18</span> <span class="n">J</span><span class="o">=</span><span class="n">I</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="mi">2</span> |
---|
461 | <a name="ln-93"></a><span class="nl"> </span> <span class="n">SCR</span><span class="o">=</span><span class="n">SCR</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">K</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
462 | <a name="ln-94"></a><span class="nl"> </span> <span class="n">SCI</span><span class="o">=</span><span class="n">SCI</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">K</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
463 | <a name="ln-95"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">2</span> |
---|
464 | <a name="ln-96"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">2</span> |
---|
465 | <a name="ln-97"></a><span class="nl">18 </span> <span class="k">CONTINUE</span> |
---|
466 | <a name="ln-98"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LLS</span><span class="o">+</span><span class="mi">1</span> |
---|
467 | <a name="ln-99"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLPS</span><span class="o">+</span><span class="mi">1</span> |
---|
468 | <a name="ln-100"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">19</span> <span class="n">J</span><span class="o">=</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="mi">2</span> |
---|
469 | <a name="ln-101"></a><span class="nl"> </span> <span class="n">ACR</span><span class="o">=</span><span class="n">ACR</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">K</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
470 | <a name="ln-102"></a><span class="nl"> </span> <span class="n">ACI</span><span class="o">=</span><span class="n">ACI</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">K</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
471 | <a name="ln-103"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">2</span> |
---|
472 | <a name="ln-104"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">2</span> |
---|
473 | <a name="ln-105"></a><span class="nl">19 </span> <span class="k">CONTINUE</span> |
---|
474 | <a name="ln-106"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
475 | <a name="ln-107"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LLS</span><span class="o">+</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">-</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span> |
---|
476 | <a name="ln-108"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLPS</span><span class="o">+</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">-</span><span class="n">I</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span> |
---|
477 | <a name="ln-109"></a><span class="nl"> </span> <span class="n">CXMS</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">SCR</span><span class="o">+</span><span class="n">ACR</span> |
---|
478 | <a name="ln-110"></a><span class="nl"> </span> <span class="n">CXMS</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">SCI</span><span class="o">+</span><span class="n">ACI</span> |
---|
479 | <a name="ln-111"></a><span class="nl"> </span> <span class="n">CXMA</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">SCR</span><span class="o">-</span><span class="n">ACR</span> |
---|
480 | <a name="ln-112"></a><span class="nl"> </span> <span class="n">CXMA</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">SCI</span><span class="o">-</span><span class="n">ACI</span> |
---|
481 | <a name="ln-113"></a><span class="nl"> 1</span> <span class="k">CONTINUE</span> |
---|
482 | <a name="ln-114"></a><span class="c">C CALL FOURTR(CXMS,FELD(IND(k)+1,L),WSAVE(:,K),MNAUF,</span> |
---|
483 | <a name="ln-115"></a><span class="c">C *MLAT(K),1)</span> |
---|
484 | <a name="ln-116"></a><span class="c">C CALL FOURTR(CXMA,FELD(MAXL-IND(k)-MLAT(K)+1,L),</span> |
---|
485 | <a name="ln-117"></a><span class="c">C *WSAVE(:,K),MNAUF,MLAT(K),1)</span> |
---|
486 | <a name="ln-118"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXMS</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">(:,</span><span class="n">K</span><span class="p">),</span><span class="n">IFAX</span><span class="p">(:,</span><span class="n">K</span><span class="p">),</span><span class="n">MNAUF</span><span class="p">,</span> |
---|
487 | <a name="ln-119"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">),</span><span class="mi">1</span><span class="p">)</span> |
---|
488 | <a name="ln-120"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">IND</span><span class="p">(</span><span class="n">k</span><span class="p">)</span><span class="o">+</span><span class="mi">1</span><span class="p">:</span><span class="n">IND</span><span class="p">(</span><span class="n">K</span><span class="p">)</span><span class="o">+</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">),</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXMS</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
489 | <a name="ln-121"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXMA</span><span class="p">,</span> |
---|
490 | <a name="ln-122"></a><span class="nl"> </span><span class="gs">*</span><span class="n">WSAVE</span><span class="p">(:,</span><span class="n">K</span><span class="p">),</span><span class="n">IFAX</span><span class="p">(:,</span><span class="n">K</span><span class="p">),</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">),</span><span class="mi">1</span><span class="p">)</span> |
---|
491 | <a name="ln-123"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="o">-</span><span class="n">IND</span><span class="p">(</span><span class="n">k</span><span class="p">)</span><span class="o">-</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">)</span><span class="o">+</span><span class="mi">1</span><span class="p">:</span><span class="n">MAXL</span><span class="o">-</span><span class="n">IND</span><span class="p">(</span><span class="n">k</span><span class="p">),</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXMA</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MLAT</span><span class="p">(</span><span class="n">K</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
492 | <a name="ln-124"></a><span class="c">C WRITE(*,*) IND+1,FELD(IND+1,L)</span> |
---|
493 | <a name="ln-125"></a><span class="nl">6 </span><span class="gs">C</span><span class="n">ONTINUE</span> |
---|
494 | <a name="ln-126"></a> |
---|
495 | <a name="ln-127"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PHSYM</span> |
---|
496 | <a name="ln-128"></a> |
---|
497 | <a name="ln-129"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PHGCUT</span><span class="p">(</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span> |
---|
498 | <a name="ln-130"></a><span class="nl"> </span><span class="gs">*</span> <span class="n">MNAUF</span><span class="p">,</span><span class="n">MMAX</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MANF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
499 | <a name="ln-131"></a> |
---|
500 | <a name="ln-132"></a><span class="c">C DIE ROUTINE FUEHRT EINE TRANSFORMATION EINER</span> |
---|
501 | <a name="ln-133"></a><span class="c">C FELDVARIABLEN VOM PHASENRAUM IN DEN PHYSIKALISCHEN</span> |
---|
502 | <a name="ln-134"></a><span class="c">C RAUM AUF KUGELKOORDINATEN DURCH. Es kann ein Teilausschnitt</span> |
---|
503 | <a name="ln-135"></a><span class="c">C Der Erde angegeben werden. Diese Routine ist langsamer als</span> |
---|
504 | <a name="ln-136"></a><span class="c">C phgrph</span> |
---|
505 | <a name="ln-137"></a><span class="c">C</span> |
---|
506 | <a name="ln-138"></a><span class="c">C CXMN = SPEKTRALKOEFFIZIENTEN IN DER REIHENFOLGE</span> |
---|
507 | <a name="ln-139"></a><span class="c">C CX00,CX01,CX11,CX02,....CXMNAUFMNAUF</span> |
---|
508 | <a name="ln-140"></a><span class="c">C FELD = FELD DER METEOROLOGISCHEN VARIABLEN</span> |
---|
509 | <a name="ln-141"></a><span class="c">C BREITE = SINUS DER GEOGRAFISCHEN BREITEN</span> |
---|
510 | <a name="ln-142"></a><span class="c">C</span> |
---|
511 | <a name="ln-143"></a><span class="c">C MNAUF ANZAHL DER FOURIERKOEFFIZIENTEN</span> |
---|
512 | <a name="ln-144"></a><span class="c">C MAUF ANZAHL DER LAENGEN UND DER FOURIERKOEFFIZIENTEN</span> |
---|
513 | <a name="ln-145"></a><span class="c">C MANF ANFANG DES LAENGENBEREICHS FUER DAS GITTER,</span> |
---|
514 | <a name="ln-146"></a><span class="c">C AUF DAS INTERPOLIERT WERDEN SOLL</span> |
---|
515 | <a name="ln-147"></a><span class="c">C MAXL ANZAHL DER FUER DAS GITTER BENUTZTEN LAENGEN</span> |
---|
516 | <a name="ln-148"></a><span class="c">C MAXB ANZAHL DER FUER DAS GITTER BENOETIGTEN BREITEN</span> |
---|
517 | <a name="ln-149"></a><span class="c">C MLEVEL ANZAHL DER LEVELS, DIE TRANSFORMIERT WERDEN</span> |
---|
518 | <a name="ln-150"></a><span class="c">C</span> |
---|
519 | <a name="ln-151"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
520 | <a name="ln-152"></a> |
---|
521 | <a name="ln-153"></a><span class="c">C FELD DER FOURIERKOEFFIZIENTEN</span> |
---|
522 | <a name="ln-154"></a> |
---|
523 | <a name="ln-155"></a><span class="c">C FELD DER LEGENDREPOLYNOME FUER EINE BREITE</span> |
---|
524 | <a name="ln-156"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MMAX</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MMAX</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="p">)</span> |
---|
525 | <a name="ln-157"></a> |
---|
526 | <a name="ln-158"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MMAX</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MMAX</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
527 | <a name="ln-159"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
528 | <a name="ln-160"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">WSAVE</span><span class="p">(</span><span class="mi">4</span><span class="o">*</span><span class="n">MAUF</span><span class="o">+</span><span class="mi">15</span><span class="p">)</span> |
---|
529 | <a name="ln-161"></a><span class="nl"> </span> <span class="kt">INTEGER</span><span class="kd">::</span> <span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span> |
---|
530 | <a name="ln-162"></a><span class="nl"> </span> |
---|
531 | <a name="ln-163"></a><span class="nl"> </span> <span class="kt">LOGICAL </span><span class="n">SYM</span> |
---|
532 | <a name="ln-164"></a> |
---|
533 | <a name="ln-165"></a><span class="c">C</span> |
---|
534 | <a name="ln-166"></a><span class="c">C write(*,*)mauf,mnauf,manf,maxl</span> |
---|
535 | <a name="ln-167"></a> |
---|
536 | <a name="ln-168"></a><span class="nl"> </span> |
---|
537 | <a name="ln-169"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">MAUF</span><span class="p">.</span><span class="n">LE</span><span class="p">.</span><span class="n">MNAUF</span><span class="p">)</span> <span class="k">WRITE</span><span class="p">(</span><span class="o">*</span><span class="p">,</span><span class="o">*</span><span class="p">)</span> <span class="s1">'TOO COARSE LONGITUDE RESOLUTION'</span> |
---|
538 | <a name="ln-170"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">((</span><span class="n">MANF</span><span class="p">.</span><span class="n">LT</span><span class="p">.</span><span class="mi">1</span><span class="p">).</span><span class="nb">OR</span><span class="p">.(</span><span class="n">MAXL</span><span class="p">.</span><span class="n">LT</span><span class="p">.</span><span class="mi">1</span><span class="p">).</span><span class="nb">OR</span><span class="p">.</span> |
---|
539 | <a name="ln-171"></a><span class="nl"> </span><span class="gs">*</span> <span class="p">(</span><span class="n">MANF</span><span class="p">.</span><span class="n">GT</span><span class="p">.</span><span class="n">MAUF</span><span class="p">).</span><span class="nb">OR</span><span class="p">.(</span><span class="n">MAXL</span><span class="p">.</span><span class="n">GT</span><span class="p">.</span><span class="n">MAUF</span><span class="p">))</span> <span class="k">THEN</span> |
---|
540 | <a name="ln-172"></a><span class="nl"> </span> <span class="k">WRITE</span><span class="p">(</span><span class="o">*</span><span class="p">,</span><span class="o">*</span><span class="p">)</span> <span class="s1">'WRONG LONGITUDE RANGE'</span><span class="p">,</span><span class="n">MANF</span><span class="p">,</span><span class="n">MAXL</span> |
---|
541 | <a name="ln-173"></a><span class="nl"> </span> <span class="k">STOP</span> |
---|
542 | <a name="ln-174"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
543 | <a name="ln-175"></a> |
---|
544 | <a name="ln-176"></a><span class="c">C Pruefe, ob Ausgabegitter symmetrisch zum Aequator ist</span> |
---|
545 | <a name="ln-177"></a><span class="c">C Wenn ja soll Symmetrie der Legendrepolynome ausgenutzt werden</span> |
---|
546 | <a name="ln-178"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">MAXB</span> <span class="p">.</span><span class="n">GT</span><span class="p">.</span> <span class="mi">4</span><span class="p">)</span> <span class="k">THEN</span> |
---|
547 | <a name="ln-179"></a><span class="nl"> </span> <span class="n">SYM</span><span class="o">=</span><span class="p">.</span><span class="n">TRUE</span><span class="p">.</span> |
---|
548 | <a name="ln-180"></a><span class="nl"> </span><span class="gs">D</span><span class="n">O</span> <span class="mi">11</span> <span class="n">J</span><span class="o">=</span><span class="mi">5</span><span class="p">,</span><span class="mi">5</span> |
---|
549 | <a name="ln-181"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="nb">ABS</span><span class="p">(</span><span class="nb">ABS</span><span class="p">(</span><span class="n">Z</span><span class="p">(</span><span class="mi">100</span><span class="p">,</span><span class="n">J</span><span class="p">))</span><span class="o">-</span><span class="nb">ABS</span><span class="p">(</span><span class="n">Z</span><span class="p">(</span><span class="mi">100</span><span class="p">,</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">J</span><span class="p">))).</span><span class="n">GT</span><span class="p">.</span><span class="err">1</span><span class="n">E</span><span class="o">-</span><span class="mi">11</span><span class="p">)</span> |
---|
550 | <a name="ln-182"></a><span class="nl"> </span><span class="gs">*</span> <span class="n">SYM</span><span class="o">=</span><span class="p">.</span><span class="n">FALSE</span><span class="p">.</span> |
---|
551 | <a name="ln-183"></a><span class="c">C WRITE(*,*) ABS(Z(100,J)),ABS(Z(100,MAXB+1-J))</span> |
---|
552 | <a name="ln-184"></a><span class="nl">11 </span><span class="gs">C</span><span class="n">ONTINUE</span> |
---|
553 | <a name="ln-185"></a><span class="nl"> </span> <span class="k">WRITE</span><span class="p">(</span><span class="o">*</span><span class="p">,</span><span class="o">*</span><span class="p">)</span> <span class="s1">'Symmetrisch: '</span><span class="p">,</span><span class="n">SYM</span> |
---|
554 | <a name="ln-186"></a><span class="nl"> </span> <span class="k">ELSE</span> |
---|
555 | <a name="ln-187"></a><span class="nl"> </span> <span class="n">SYM</span><span class="o">=</span><span class="p">.</span><span class="n">FALSE</span><span class="p">.</span> |
---|
556 | <a name="ln-188"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
557 | <a name="ln-189"></a> |
---|
558 | <a name="ln-190"></a> |
---|
559 | <a name="ln-191"></a><span class="nl"> IF(</span><span class="gs">S</span><span class="n">YM</span><span class="p">)</span> <span class="k">THEN</span> |
---|
560 | <a name="ln-192"></a><span class="c">!$OMP PARALLEL DO </span> |
---|
561 | <a name="ln-193"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">J</span><span class="o">=</span><span class="mi">1</span><span class="p">,(</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">/</span><span class="mi">2</span> |
---|
562 | <a name="ln-194"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">PHSYMCUT</span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span> |
---|
563 | <a name="ln-195"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">,</span><span class="n">MANF</span><span class="p">)</span> |
---|
564 | <a name="ln-196"></a> |
---|
565 | <a name="ln-197"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
566 | <a name="ln-198"></a><span class="c">!$OMP END PARALLEL DO</span> |
---|
567 | <a name="ln-199"></a><span class="nl"> ELS</span><span class="gs">E</span> |
---|
568 | <a name="ln-200"></a><span class="c">!$OMP PARALLEL DO </span> |
---|
569 | <a name="ln-201"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">J</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXB</span> |
---|
570 | <a name="ln-202"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">PHGPNS</span><span class="p">(</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span> |
---|
571 | <a name="ln-203"></a><span class="nl"> </span><span class="gs">*</span><span class="n">J</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MANF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
572 | <a name="ln-204"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
573 | <a name="ln-205"></a><span class="c">!$OMP END PARALLEL DO</span> |
---|
574 | <a name="ln-206"></a> |
---|
575 | <a name="ln-207"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
576 | <a name="ln-208"></a> |
---|
577 | <a name="ln-209"></a> |
---|
578 | <a name="ln-210"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
579 | <a name="ln-211"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PHGCUT</span> |
---|
580 | <a name="ln-212"></a> |
---|
581 | <a name="ln-213"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PHSYMCUT</span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span> |
---|
582 | <a name="ln-214"></a><span class="nl"> </span><span class="gs">*</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">,</span><span class="n">MANF</span><span class="p">)</span> |
---|
583 | <a name="ln-215"></a> |
---|
584 | <a name="ln-216"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
585 | <a name="ln-217"></a> |
---|
586 | <a name="ln-218"></a><span class="c">C FELD DER FOURIERKOEFFIZIENTEN</span> |
---|
587 | <a name="ln-219"></a> |
---|
588 | <a name="ln-220"></a><span class="nl"> </span> <span class="kt">REAL</span> <span class="kd">::</span> <span class="n">CXM</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">),</span><span class="n">CXMA</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
589 | <a name="ln-221"></a> |
---|
590 | <a name="ln-222"></a> |
---|
591 | <a name="ln-223"></a><span class="c">C FELD DER LEGENDREPOLYNOME FUER EINE BREITE</span> |
---|
592 | <a name="ln-224"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="p">)</span> |
---|
593 | <a name="ln-225"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">SCR</span><span class="p">,</span><span class="n">SCI</span><span class="p">,</span><span class="n">ACR</span><span class="p">,</span><span class="n">ACI</span> |
---|
594 | <a name="ln-226"></a> |
---|
595 | <a name="ln-227"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
596 | <a name="ln-228"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
597 | <a name="ln-229"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">WSAVE</span><span class="p">(</span><span class="mi">4</span><span class="o">*</span><span class="n">MAUF</span><span class="o">+</span><span class="mi">15</span><span class="p">)</span> |
---|
598 | <a name="ln-230"></a><span class="nl"> </span> <span class="kt">INTEGER</span> <span class="kd">::</span> <span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span> |
---|
599 | <a name="ln-231"></a> |
---|
600 | <a name="ln-232"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">16</span> <span class="n">L</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span> |
---|
601 | <a name="ln-233"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="mi">0</span> |
---|
602 | <a name="ln-234"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="mi">0</span> |
---|
603 | <a name="ln-235"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">17</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
604 | <a name="ln-236"></a><span class="nl"> </span> <span class="n">SCR</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
605 | <a name="ln-237"></a><span class="nl"> </span> <span class="n">SCI</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
606 | <a name="ln-238"></a><span class="nl"> </span> <span class="n">ACR</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
607 | <a name="ln-239"></a><span class="nl"> </span> <span class="n">ACI</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
608 | <a name="ln-240"></a><span class="nl"> </span> <span class="n">LLS</span><span class="o">=</span><span class="n">LL</span> |
---|
609 | <a name="ln-241"></a><span class="nl"> </span> <span class="n">LLPS</span><span class="o">=</span><span class="n">LLP</span> |
---|
610 | <a name="ln-242"></a><span class="c">C Innerste Schleife aufgespalten um if-Abfrage zu sparen</span> |
---|
611 | <a name="ln-243"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">18</span> <span class="n">K</span><span class="o">=</span><span class="n">I</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="mi">2</span> |
---|
612 | <a name="ln-244"></a><span class="nl"> </span> <span class="n">SCR</span><span class="o">=</span><span class="n">SCR</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
613 | <a name="ln-245"></a><span class="nl"> </span> <span class="n">SCI</span><span class="o">=</span><span class="n">SCI</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
614 | <a name="ln-246"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">2</span> |
---|
615 | <a name="ln-247"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">2</span> |
---|
616 | <a name="ln-248"></a><span class="nl">18 </span> <span class="k">CONTINUE</span> |
---|
617 | <a name="ln-249"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LLS</span><span class="o">+</span><span class="mi">1</span> |
---|
618 | <a name="ln-250"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLPS</span><span class="o">+</span><span class="mi">1</span> |
---|
619 | <a name="ln-251"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">19</span> <span class="n">K</span><span class="o">=</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="mi">2</span> |
---|
620 | <a name="ln-252"></a><span class="nl"> </span> <span class="n">ACR</span><span class="o">=</span><span class="n">ACR</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
621 | <a name="ln-253"></a><span class="nl"> </span> <span class="n">ACI</span><span class="o">=</span><span class="n">ACI</span> <span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">,</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">L</span><span class="p">)</span> |
---|
622 | <a name="ln-254"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">2</span> |
---|
623 | <a name="ln-255"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">2</span> |
---|
624 | <a name="ln-256"></a><span class="nl">19 </span> <span class="k">CONTINUE</span> |
---|
625 | <a name="ln-257"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LLS</span><span class="o">+</span><span class="n">MNAUF</span><span class="o">-</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span> |
---|
626 | <a name="ln-258"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLPS</span><span class="o">+</span><span class="n">MNAUF</span><span class="o">-</span><span class="n">I</span><span class="o">+</span><span class="mi">3</span> |
---|
627 | <a name="ln-259"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">SCR</span><span class="o">+</span><span class="n">ACR</span> |
---|
628 | <a name="ln-260"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">SCI</span><span class="o">+</span><span class="n">ACI</span> |
---|
629 | <a name="ln-261"></a><span class="nl"> </span> <span class="n">CXMA</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">SCR</span><span class="o">-</span><span class="n">ACR</span> |
---|
630 | <a name="ln-262"></a><span class="nl"> </span> <span class="n">CXMA</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">SCI</span><span class="o">-</span><span class="n">ACI</span> |
---|
631 | <a name="ln-263"></a><span class="nl">17 </span> <span class="k">CONTINUE</span> |
---|
632 | <a name="ln-264"></a><span class="nl"> </span> |
---|
633 | <a name="ln-265"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXM</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span> |
---|
634 | <a name="ln-266"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">26</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MAXL</span><span class="o">-</span><span class="mi">1</span> |
---|
635 | <a name="ln-267"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="p">.</span><span class="n">LE</span><span class="p">.</span><span class="n">MAUF</span><span class="p">)</span> <span class="k">THEN</span> |
---|
636 | <a name="ln-268"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXM</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
637 | <a name="ln-269"></a><span class="nl"> </span> <span class="k">ELSE</span> |
---|
638 | <a name="ln-270"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXM</span><span class="p">(</span><span class="n">MANF</span><span class="o">-</span><span class="n">MAUF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
639 | <a name="ln-271"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
640 | <a name="ln-272"></a><span class="nl">26 </span> <span class="k">CONTINUE</span> |
---|
641 | <a name="ln-273"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXMA</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span> |
---|
642 | <a name="ln-274"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">36</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MAXL</span><span class="o">-</span><span class="mi">1</span> |
---|
643 | <a name="ln-275"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="p">.</span><span class="n">LE</span><span class="p">.</span><span class="n">MAUF</span><span class="p">)</span> <span class="k">THEN</span> |
---|
644 | <a name="ln-276"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXMA</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
645 | <a name="ln-277"></a><span class="nl"> </span> <span class="k">ELSE</span> |
---|
646 | <a name="ln-278"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXB</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXMA</span><span class="p">(</span><span class="n">MANF</span><span class="o">-</span><span class="n">MAUF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
647 | <a name="ln-279"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
648 | <a name="ln-280"></a><span class="nl">36 </span> <span class="k">CONTINUE</span> |
---|
649 | <a name="ln-281"></a><span class="nl">16 </span><span class="gs"> </span><span class="k">CONTINUE</span> |
---|
650 | <a name="ln-282"></a> |
---|
651 | <a name="ln-283"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PHSYMCUT</span> |
---|
652 | <a name="ln-284"></a> |
---|
653 | <a name="ln-285"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PHGPNS</span><span class="p">(</span><span class="n">CXMN</span><span class="p">,</span><span class="n">FELD</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span> |
---|
654 | <a name="ln-286"></a><span class="nl"> </span><span class="gs">*</span><span class="n">J</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MANF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
655 | <a name="ln-287"></a> |
---|
656 | <a name="ln-288"></a><span class="nl"> </span> <span class="k">IMPLICIT NONE</span> |
---|
657 | <a name="ln-289"></a><span class="nl"> </span> <span class="kt">INTEGER</span><span class="p">,</span><span class="k">intent</span><span class="p">(</span><span class="n">in</span><span class="p">)</span> <span class="kd">::</span> <span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="n">MANF</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span> |
---|
658 | <a name="ln-290"></a><span class="nl"> </span> <span class="kt">REAL</span> <span class="kd">::</span> <span class="n">CXM</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
659 | <a name="ln-291"></a><span class="nl"> </span> <span class="kt">REAL</span><span class="p">,</span><span class="k">intent</span><span class="p">(</span><span class="n">in</span><span class="p">)</span> <span class="kd">::</span> <span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">,</span><span class="n">MAXB</span><span class="p">)</span> |
---|
660 | <a name="ln-292"></a> |
---|
661 | <a name="ln-293"></a><span class="nl"> </span> <span class="kt">REAL</span><span class="p">,</span><span class="k">intent</span><span class="p">(</span><span class="n">in</span><span class="p">)</span> <span class="kd">::</span> <span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
662 | <a name="ln-294"></a> |
---|
663 | <a name="ln-295"></a><span class="nl"> </span> <span class="kt">REAL</span><span class="p">,</span><span class="k">intent</span><span class="p">(</span><span class="n">in</span><span class="p">)</span> <span class="kd">::</span> <span class="n">WSAVE</span><span class="p">(</span><span class="mi">4</span><span class="o">*</span><span class="n">MAUF</span><span class="o">+</span><span class="mi">15</span><span class="p">)</span> |
---|
664 | <a name="ln-296"></a> |
---|
665 | <a name="ln-297"></a><span class="nl"> </span> <span class="kt">REAL</span> <span class="kd">::</span> <span class="n">FELD</span><span class="p">(</span><span class="n">MAXL</span><span class="p">,</span><span class="n">MAXB</span><span class="p">,</span><span class="n">MLEVEL</span><span class="p">)</span> |
---|
666 | <a name="ln-298"></a><span class="nl"> </span> <span class="kt">INTEGER</span> <span class="kd">::</span> <span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span> |
---|
667 | <a name="ln-299"></a> |
---|
668 | <a name="ln-300"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">I</span><span class="p">,</span><span class="n">L</span> |
---|
669 | <a name="ln-301"></a> |
---|
670 | <a name="ln-302"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">L</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">MLEVEL</span> |
---|
671 | <a name="ln-303"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">LEGTR</span><span class="p">(</span><span class="n">CXMN</span><span class="p">(:,</span><span class="n">L</span><span class="p">),</span><span class="n">CXM</span><span class="p">,</span><span class="n">Z</span><span class="p">(:,</span><span class="n">J</span><span class="p">),</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">)</span> |
---|
672 | <a name="ln-304"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXM</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span> |
---|
673 | <a name="ln-305"></a><span class="nl"> </span> |
---|
674 | <a name="ln-306"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MAXL</span><span class="o">-</span><span class="mi">1</span> |
---|
675 | <a name="ln-307"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="p">.</span><span class="n">LE</span><span class="p">.</span><span class="n">MAUF</span><span class="p">)</span> <span class="k">THEN</span> |
---|
676 | <a name="ln-308"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXM</span><span class="p">(</span><span class="n">MANF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
677 | <a name="ln-309"></a><span class="nl"> </span> <span class="k">ELSE</span> |
---|
678 | <a name="ln-310"></a><span class="nl"> </span> <span class="n">FELD</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">L</span><span class="p">)</span><span class="o">=</span><span class="n">CXM</span><span class="p">(</span><span class="n">MANF</span><span class="o">-</span><span class="n">MAUF</span><span class="o">+</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
679 | <a name="ln-311"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
680 | <a name="ln-312"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
681 | <a name="ln-313"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
682 | <a name="ln-314"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PHGPNS</span> |
---|
683 | <a name="ln-315"></a><span class="c">C</span> |
---|
684 | <a name="ln-316"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">LEGTR</span><span class="p">(</span><span class="n">CXMN</span><span class="p">,</span><span class="n">CXM</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">)</span> |
---|
685 | <a name="ln-317"></a><span class="nl"> </span> <span class="k">IMPLICIT NONE</span> |
---|
686 | <a name="ln-318"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAUF</span><span class="p">,</span><span class="n">LL</span><span class="p">,</span><span class="n">LLP</span><span class="p">,</span><span class="n">I</span><span class="p">,</span><span class="n">J</span> |
---|
687 | <a name="ln-319"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">CXM</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
688 | <a name="ln-320"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">CXMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
689 | <a name="ln-321"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
690 | <a name="ln-322"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">CI</span><span class="p">,</span><span class="n">CR</span> |
---|
691 | <a name="ln-323"></a><span class="c">C</span> |
---|
692 | <a name="ln-324"></a><span class="c">C DIESE ROUTINE BERECHNET DIE FOURIERKOEFFIZIENTEN CXM</span> |
---|
693 | <a name="ln-325"></a><span class="c">C</span> |
---|
694 | <a name="ln-326"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="mi">0</span> |
---|
695 | <a name="ln-327"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="mi">0</span> |
---|
696 | <a name="ln-328"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">1</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
697 | <a name="ln-329"></a><span class="nl"> </span> <span class="n">CR</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
698 | <a name="ln-330"></a><span class="nl"> </span> <span class="n">CI</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
699 | <a name="ln-331"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">2</span> <span class="n">J</span><span class="o">=</span><span class="n">I</span><span class="p">,</span><span class="n">MNAUF</span> |
---|
700 | <a name="ln-332"></a><span class="nl"> </span> <span class="n">CR</span><span class="o">=</span><span class="n">CR</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="p">)</span> |
---|
701 | <a name="ln-333"></a><span class="nl"> </span> <span class="n">CI</span><span class="o">=</span><span class="n">CI</span><span class="o">+</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="p">)</span><span class="o">*</span><span class="n">CXMN</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span> |
---|
702 | <a name="ln-334"></a><span class="nl"> </span> <span class="n">LL</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">1</span> |
---|
703 | <a name="ln-335"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">1</span> |
---|
704 | <a name="ln-336"></a><span class="nl"> 2</span> <span class="k">CONTINUE</span> |
---|
705 | <a name="ln-337"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">2</span> |
---|
706 | <a name="ln-338"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">CR</span> |
---|
707 | <a name="ln-339"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">CI</span> |
---|
708 | <a name="ln-340"></a><span class="nl"> 1</span> <span class="k">CONTINUE</span> |
---|
709 | <a name="ln-341"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
710 | <a name="ln-342"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">LEGTR</span> |
---|
711 | <a name="ln-343"></a><span class="c">C</span> |
---|
712 | <a name="ln-344"></a><span class="c">C </span> |
---|
713 | <a name="ln-345"></a><span class="c">C</span> |
---|
714 | <a name="ln-346"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">RFOURTR</span><span class="p">(</span><span class="n">CXM</span><span class="p">,</span><span class="n">TRIGS</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="nb">ISIGN</span><span class="p">)</span> |
---|
715 | <a name="ln-347"></a><span class="c">C BERECHNET DIE FOURIERSUMME MIT EINEM FFT-ALGORITHMUS</span> |
---|
716 | <a name="ln-348"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
717 | <a name="ln-349"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">CXM</span><span class="p">(</span><span class="mi">0</span><span class="p">:</span><span class="n">MAXAUF</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
718 | <a name="ln-350"></a><span class="nl"> </span> <span class="kt">REAL</span> <span class="kd">::</span> <span class="n">WSAVE</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">MAXL</span><span class="p">),</span><span class="n">TRIGS</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">MAXL</span><span class="p">)</span> |
---|
719 | <a name="ln-351"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">IFAX</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span> |
---|
720 | <a name="ln-352"></a> |
---|
721 | <a name="ln-353"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">I</span><span class="o">=</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXL</span><span class="o">-</span><span class="mi">1</span> |
---|
722 | <a name="ln-354"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="mf">0.0</span> |
---|
723 | <a name="ln-355"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="mf">0.0</span> |
---|
724 | <a name="ln-356"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
725 | <a name="ln-357"></a><span class="nl"> </span> <span class="k">CALL </span><span class="n">FFT99</span><span class="p">(</span><span class="n">CXM</span><span class="p">,</span><span class="n">WSAVE</span><span class="p">,</span><span class="n">TRIGS</span><span class="p">,</span><span class="n">IFAX</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="n">MAXL</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span> |
---|
726 | <a name="ln-358"></a><span class="nl"> </span> <span class="k">DO </span><span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MAXL</span><span class="o">-</span><span class="mi">1</span> |
---|
727 | <a name="ln-359"></a><span class="nl"> </span> <span class="n">CXM</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">CXM</span><span class="p">(</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span> |
---|
728 | <a name="ln-360"></a><span class="nl"> </span> <span class="n">ENDDO</span> |
---|
729 | <a name="ln-361"></a> |
---|
730 | <a name="ln-362"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
731 | <a name="ln-363"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">RFOURTR</span> |
---|
732 | <a name="ln-364"></a><span class="c">C </span> |
---|
733 | <a name="ln-365"></a><span class="c">C </span> |
---|
734 | <a name="ln-366"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">GAULEG</span><span class="p">(</span><span class="n">X1</span><span class="p">,</span><span class="n">X2</span><span class="p">,</span><span class="n">X</span><span class="p">,</span><span class="n">W</span><span class="p">,</span><span class="n">N</span><span class="p">)</span> |
---|
735 | <a name="ln-367"></a><span class="c">C BERECHNET DIE GAUSS+SCHEN BREITEN</span> |
---|
736 | <a name="ln-368"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
737 | <a name="ln-369"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">X</span><span class="p">(</span><span class="n">N</span><span class="p">),</span><span class="n">W</span><span class="p">(</span><span class="n">N</span><span class="p">)</span> |
---|
738 | <a name="ln-370"></a><span class="nl"> </span> <span class="k">PARAMETER</span> <span class="p">(</span><span class="n">EPS</span><span class="o">=</span><span class="mf">3.D-14</span><span class="p">)</span> |
---|
739 | <a name="ln-371"></a><span class="nl"> </span> <span class="n">M</span><span class="o">=</span><span class="p">(</span><span class="n">N</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">/</span><span class="mi">2</span> |
---|
740 | <a name="ln-372"></a><span class="nl"> </span> <span class="n">XM</span><span class="o">=</span><span class="mf">0.5D0</span><span class="o">*</span><span class="p">(</span><span class="n">X2</span><span class="o">+</span><span class="n">X1</span><span class="p">)</span> |
---|
741 | <a name="ln-373"></a><span class="nl"> </span> <span class="n">XL</span><span class="o">=</span><span class="mf">0.5D0</span><span class="o">*</span><span class="p">(</span><span class="n">X2</span><span class="o">-</span><span class="n">X1</span><span class="p">)</span> |
---|
742 | <a name="ln-374"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">12</span> <span class="n">I</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">M</span> |
---|
743 | <a name="ln-375"></a><span class="nl"> </span> <span class="n">Z</span><span class="o">=</span><span class="nb">DCOS</span><span class="p">(</span><span class="mf">3.141592654D0</span><span class="o">*</span><span class="p">(</span><span class="n">I</span><span class="o">-</span><span class="p">.</span><span class="mi">25</span><span class="n">D0</span><span class="p">)</span><span class="o">/</span><span class="p">(</span><span class="n">N</span><span class="o">+</span><span class="p">.</span><span class="mi">5</span><span class="n">D0</span><span class="p">))</span> |
---|
744 | <a name="ln-376"></a><span class="nl">1 </span> <span class="k">CONTINUE</span> |
---|
745 | <a name="ln-377"></a><span class="nl"> </span> <span class="n">P1</span><span class="o">=</span><span class="mf">1.D0</span> |
---|
746 | <a name="ln-378"></a><span class="nl"> </span> <span class="n">P2</span><span class="o">=</span><span class="mf">0.D0</span> |
---|
747 | <a name="ln-379"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">11</span> <span class="n">J</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">N</span> |
---|
748 | <a name="ln-380"></a><span class="nl"> </span> <span class="n">P3</span><span class="o">=</span><span class="n">P2</span> |
---|
749 | <a name="ln-381"></a><span class="nl"> </span> <span class="n">P2</span><span class="o">=</span><span class="n">P1</span> |
---|
750 | <a name="ln-382"></a><span class="nl"> </span> <span class="n">P1</span><span class="o">=</span><span class="p">((</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">J</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">*</span><span class="n">Z</span><span class="o">*</span><span class="n">P2</span><span class="o">-</span><span class="p">(</span><span class="n">J</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">*</span><span class="n">P3</span><span class="p">)</span><span class="o">/</span><span class="n">J</span> |
---|
751 | <a name="ln-383"></a><span class="nl">11 </span> <span class="k">CONTINUE</span> |
---|
752 | <a name="ln-384"></a><span class="nl"> </span> <span class="n">PP</span><span class="o">=</span><span class="n">N</span><span class="o">*</span><span class="p">(</span><span class="n">Z</span><span class="o">*</span><span class="n">P1</span><span class="o">-</span><span class="n">P2</span><span class="p">)</span><span class="o">/</span><span class="p">(</span><span class="n">Z</span><span class="o">*</span><span class="n">Z</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span> |
---|
753 | <a name="ln-385"></a><span class="nl"> </span> <span class="n">Z1</span><span class="o">=</span><span class="n">Z</span> |
---|
754 | <a name="ln-386"></a><span class="nl"> </span> <span class="n">Z</span><span class="o">=</span><span class="n">Z1</span><span class="o">-</span><span class="n">P1</span><span class="o">/</span><span class="n">PP</span> |
---|
755 | <a name="ln-387"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="nb">ABS</span><span class="p">(</span><span class="n">Z</span><span class="o">-</span><span class="n">Z1</span><span class="p">).</span><span class="n">GT</span><span class="p">.</span><span class="n">EPS</span><span class="p">)</span><span class="n">GO</span> <span class="n">TO</span> <span class="mi">1</span> |
---|
756 | <a name="ln-388"></a><span class="nl"> </span> <span class="n">X</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">XM</span><span class="o">-</span><span class="n">XL</span><span class="o">*</span><span class="n">Z</span> |
---|
757 | <a name="ln-389"></a><span class="nl"> </span> <span class="n">X</span><span class="p">(</span><span class="n">N</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">XM</span><span class="o">+</span><span class="n">XL</span><span class="o">*</span><span class="n">Z</span> |
---|
758 | <a name="ln-390"></a><span class="nl"> </span> <span class="n">W</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">XL</span><span class="o">/</span><span class="p">((</span><span class="mf">1.D0</span><span class="o">-</span><span class="n">Z</span><span class="o">*</span><span class="n">Z</span><span class="p">)</span><span class="o">*</span><span class="n">PP</span><span class="o">*</span><span class="n">PP</span><span class="p">)</span> |
---|
759 | <a name="ln-391"></a><span class="nl"> </span> <span class="n">W</span><span class="p">(</span><span class="n">N</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">W</span><span class="p">(</span><span class="n">I</span><span class="p">)</span> |
---|
760 | <a name="ln-392"></a><span class="nl">12 </span> <span class="k">CONTINUE</span> |
---|
761 | <a name="ln-393"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
762 | <a name="ln-394"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">GAULEG</span> |
---|
763 | <a name="ln-395"></a><span class="c">C</span> |
---|
764 | <a name="ln-396"></a><span class="c">C</span> |
---|
765 | <a name="ln-397"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">PLGNFA</span><span class="p">(</span><span class="n">LL</span><span class="p">,</span><span class="n">X</span><span class="p">,</span><span class="n">Z</span><span class="p">)</span> |
---|
766 | <a name="ln-398"></a><span class="c">C</span> |
---|
767 | <a name="ln-399"></a><span class="c">C PLGNDN BERECHNET ALLE NORMIERTEN ASSOZIIERTEN</span> |
---|
768 | <a name="ln-400"></a><span class="c">C LEGENDREFUNKTIONEN VON P00(X) BIS PLL(X)</span> |
---|
769 | <a name="ln-401"></a><span class="c">C UND SCHREIBT SIE IN DAS FELD Z</span> |
---|
770 | <a name="ln-402"></a><span class="c">C Die Polynome sind wie im ECMWF indiziert, d.h.</span> |
---|
771 | <a name="ln-403"></a><span class="c">C P00,P10,P11,P20,P21,P22,...</span> |
---|
772 | <a name="ln-404"></a><span class="c">C Ansonsten ist die Routine analog zu PLGNDN</span> |
---|
773 | <a name="ln-405"></a><span class="c">C X IST DER COSINUS DES ZENITWINKELS ODER</span> |
---|
774 | <a name="ln-406"></a><span class="c">C DER SINUS DER GEOGRAFISCHEN BREITE</span> |
---|
775 | <a name="ln-407"></a><span class="c">C</span> |
---|
776 | <a name="ln-408"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
777 | <a name="ln-409"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">LL</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">LL</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
778 | <a name="ln-410"></a><span class="c">C</span> |
---|
779 | <a name="ln-411"></a><span class="nl"> </span> <span class="n">L</span><span class="o">=</span><span class="n">LL</span><span class="o">+</span><span class="mi">2</span> |
---|
780 | <a name="ln-412"></a><span class="nl"> </span> <span class="n">I</span><span class="o">=</span><span class="mi">1</span> |
---|
781 | <a name="ln-413"></a><span class="nl"> </span> <span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">)</span><span class="o">=</span><span class="mf">1.D0</span> |
---|
782 | <a name="ln-414"></a><span class="nl"> </span> <span class="n">FACT</span><span class="o">=</span><span class="mf">1.D0</span> |
---|
783 | <a name="ln-415"></a><span class="nl"> </span> <span class="n">POT</span><span class="o">=</span><span class="mf">1.D0</span> |
---|
784 | <a name="ln-416"></a><span class="nl"> </span> <span class="n">SOMX2</span><span class="o">=</span><span class="nb">DSQRT</span><span class="p">(</span><span class="mf">1.D0</span><span class="o">-</span><span class="n">X</span><span class="o">*</span><span class="n">X</span><span class="p">)</span> |
---|
785 | <a name="ln-417"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">14</span> <span class="n">J</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">L</span> |
---|
786 | <a name="ln-418"></a><span class="nl"> </span> <span class="n">DJ</span><span class="o">=</span><span class="nb">DBLE</span><span class="p">(</span><span class="n">J</span><span class="p">)</span> |
---|
787 | <a name="ln-419"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">J</span><span class="p">.</span><span class="n">GT</span><span class="p">.</span><span class="mi">0</span><span class="p">)</span> <span class="k">THEN</span> |
---|
788 | <a name="ln-420"></a><span class="nl"> </span> <span class="n">FACT</span><span class="o">=</span><span class="n">FACT</span><span class="o">*</span><span class="p">(</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DJ</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">/</span><span class="p">(</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DJ</span><span class="p">)</span> |
---|
789 | <a name="ln-421"></a><span class="nl"> </span> <span class="n">POT</span><span class="o">=</span><span class="n">POT</span><span class="o">*</span><span class="n">SOMX2</span> |
---|
790 | <a name="ln-422"></a><span class="nl"> </span> <span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="nb">DSQRT</span><span class="p">((</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DJ</span><span class="o">+</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">*</span><span class="n">FACT</span><span class="p">)</span><span class="o">*</span><span class="n">POT</span> |
---|
791 | <a name="ln-423"></a><span class="nl"> </span> <span class="n">I</span><span class="o">=</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span> |
---|
792 | <a name="ln-424"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
793 | <a name="ln-425"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">J</span><span class="p">.</span><span class="n">LT</span><span class="p">.</span><span class="n">L</span><span class="p">)</span> <span class="k">THEN</span> |
---|
794 | <a name="ln-426"></a><span class="nl"> </span> <span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">X</span><span class="o">*</span> |
---|
795 | <a name="ln-427"></a><span class="nl"> </span><span class="gs">*</span><span class="nb">DSQRT</span><span class="p">((</span><span class="mf">4.D0</span><span class="o">*</span><span class="n">DJ</span><span class="o">*</span><span class="n">DJ</span><span class="o">+</span><span class="mf">8.D0</span><span class="o">*</span><span class="n">DJ</span><span class="o">+</span><span class="mf">3.D0</span><span class="p">)</span><span class="o">/</span><span class="p">(</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DJ</span><span class="o">+</span><span class="mf">1.D0</span><span class="p">))</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
796 | <a name="ln-428"></a><span class="nl"> </span> <span class="n">I</span><span class="o">=</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span> |
---|
797 | <a name="ln-429"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
798 | <a name="ln-430"></a><span class="nl"> </span> <span class="n">DK</span><span class="o">=</span><span class="n">DJ</span><span class="o">+</span><span class="mf">2.D0</span> |
---|
799 | <a name="ln-431"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">14</span> <span class="n">K</span><span class="o">=</span><span class="n">J</span><span class="o">+</span><span class="mi">2</span><span class="p">,</span><span class="n">L</span> |
---|
800 | <a name="ln-432"></a><span class="nl"> </span> <span class="n">DDK</span><span class="o">=</span><span class="p">(</span><span class="n">DK</span><span class="o">*</span><span class="n">DK</span><span class="o">-</span><span class="n">DJ</span><span class="o">*</span><span class="n">DJ</span><span class="p">)</span> |
---|
801 | <a name="ln-433"></a><span class="nl"> </span> <span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="p">)</span><span class="o">=</span><span class="n">X</span><span class="o">*</span><span class="nb">DSQRT</span><span class="p">((</span><span class="mf">4.D0</span><span class="o">*</span><span class="n">DK</span><span class="o">*</span><span class="n">DK</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">/</span><span class="n">DDK</span><span class="p">)</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span> |
---|
802 | <a name="ln-434"></a><span class="nl"> </span><span class="gs">*</span> <span class="nb">DSQRT</span><span class="p">(((</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DK</span><span class="o">+</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">DK</span><span class="o">-</span><span class="n">DJ</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">DK</span><span class="o">+</span><span class="n">DJ</span><span class="o">-</span><span class="mf">1.D0</span><span class="p">))</span><span class="o">/</span> |
---|
803 | <a name="ln-435"></a><span class="nl"> </span><span class="gs">*</span> <span class="p">((</span><span class="mf">2.D0</span><span class="o">*</span><span class="n">DK</span><span class="o">-</span><span class="mf">3.D0</span><span class="p">)</span><span class="o">*</span><span class="n">DDK</span><span class="p">))</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">I</span><span class="o">-</span><span class="mi">2</span><span class="p">)</span> |
---|
804 | <a name="ln-436"></a><span class="nl"> </span> <span class="n">DK</span><span class="o">=</span><span class="n">DK</span><span class="o">+</span><span class="mf">1.D0</span> |
---|
805 | <a name="ln-437"></a><span class="nl"> </span> <span class="n">I</span><span class="o">=</span><span class="n">I</span><span class="o">+</span><span class="mi">1</span> |
---|
806 | <a name="ln-438"></a><span class="nl">14 </span> <span class="k">CONTINUE</span> |
---|
807 | <a name="ln-439"></a><span class="nl"> </span> <span class="k">RETURN</span> |
---|
808 | <a name="ln-440"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">PLGNFA</span> |
---|
809 | <a name="ln-441"></a> |
---|
810 | <a name="ln-442"></a> |
---|
811 | <a name="ln-443"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">DPLGND</span><span class="p">(</span><span class="n">MNAUF</span><span class="p">,</span><span class="n">Z</span><span class="p">,</span><span class="n">DZ</span><span class="p">)</span> |
---|
812 | <a name="ln-444"></a><span class="c">C</span> |
---|
813 | <a name="ln-445"></a><span class="c">C DPLGND BERECHNET DIE ABLEITUNG DER NORMIERTEN ASSOZIIERTEN</span> |
---|
814 | <a name="ln-446"></a><span class="c">C LEGENDREFUNKTIONEN VON P00(X) BIS PLL(X)</span> |
---|
815 | <a name="ln-447"></a><span class="c">C UND SCHREIBT SIE IN DAS FELD DZ</span> |
---|
816 | <a name="ln-448"></a><span class="c">C DIE REIHENFOLGE IST</span> |
---|
817 | <a name="ln-449"></a><span class="c">C P00(X),P01(X),P11(X),P02(X),P12(X),P22(X),..PLL(X)</span> |
---|
818 | <a name="ln-450"></a><span class="c">C</span> |
---|
819 | <a name="ln-451"></a><span class="nl"> </span> <span class="k">IMPLICIT </span><span class="kt">REAL</span> <span class="p">(</span><span class="n">A</span><span class="o">-</span><span class="n">H</span><span class="p">,</span><span class="n">O</span><span class="o">-</span><span class="n">Z</span><span class="p">)</span> |
---|
820 | <a name="ln-452"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">4</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
821 | <a name="ln-453"></a><span class="nl"> </span> <span class="k">DIMENSION </span><span class="n">DZ</span><span class="p">(</span><span class="mi">0</span><span class="p">:((</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">3</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span> |
---|
822 | <a name="ln-454"></a><span class="c">C</span> |
---|
823 | <a name="ln-455"></a><span class="nl"> I</span><span class="gs">F</span><span class="p">(</span><span class="n">Z</span><span class="p">(</span><span class="mi">0</span><span class="p">).</span><span class="n">NE</span><span class="p">.</span><span class="mf">1.D0</span><span class="p">)</span> <span class="k">THEN</span> |
---|
824 | <a name="ln-456"></a><span class="nl"> WR</span><span class="gs">I</span><span class="n">TE</span><span class="p">(</span><span class="o">*</span><span class="p">,</span><span class="o">*</span><span class="p">)</span> <span class="s1">'DPLGND: Z(0) must be 1.0'</span> |
---|
825 | <a name="ln-457"></a><span class="nl"> </span><span class="gs">S</span><span class="n">TOP</span> |
---|
826 | <a name="ln-458"></a><span class="nl"> E</span><span class="gs">N</span><span class="n">DIF</span> |
---|
827 | <a name="ln-459"></a> |
---|
828 | <a name="ln-460"></a><span class="nl"> L</span><span class="gs">L</span><span class="n">P</span><span class="o">=</span><span class="mi">0</span> |
---|
829 | <a name="ln-461"></a><span class="nl"> L</span><span class="gs">L</span><span class="n">H</span><span class="o">=</span><span class="mi">0</span> |
---|
830 | <a name="ln-462"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">1</span> <span class="n">I</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span> |
---|
831 | <a name="ln-463"></a><span class="nl"> </span> <span class="k">DO </span><span class="mi">2</span> <span class="n">J</span><span class="o">=</span><span class="n">I</span><span class="p">,</span><span class="n">MNAUF</span><span class="o">+</span><span class="mi">1</span> |
---|
832 | <a name="ln-464"></a><span class="nl"> </span> <span class="k">IF</span><span class="p">(</span><span class="n">I</span><span class="p">.</span><span class="n">EQ</span><span class="p">.</span><span class="n">J</span><span class="p">)</span> <span class="k">THEN</span> |
---|
833 | <a name="ln-465"></a><span class="nl"> </span> <span class="n">WURZELA</span><span class="o">=</span> |
---|
834 | <a name="ln-466"></a><span class="nl"> </span><span class="gs">*</span><span class="nb">DSQRT</span><span class="p">(</span><span class="nb">DBLE</span><span class="p">((</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span><span class="n">I</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">/</span><span class="nb">DBLE</span><span class="p">(</span><span class="mi">4</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">))</span> |
---|
835 | <a name="ln-467"></a><span class="nl"> </span> <span class="n">DZ</span><span class="p">(</span><span class="n">LLH</span><span class="p">)</span><span class="o">=</span><span class="nb">DBLE</span><span class="p">(</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">WURZELA</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span> |
---|
836 | <a name="ln-468"></a><span class="nl"> </span> <span class="k">ELSE</span> |
---|
837 | <a name="ln-469"></a><span class="nl"> </span> <span class="n">WURZELB</span><span class="o">=</span> |
---|
838 | <a name="ln-470"></a><span class="nl"> </span><span class="gs">*</span><span class="nb">DSQRT</span><span class="p">(</span><span class="nb">DBLE</span><span class="p">((</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span><span class="n">I</span><span class="o">*</span><span class="n">I</span><span class="p">)</span><span class="o">/</span><span class="nb">DBLE</span><span class="p">(</span><span class="mi">4</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">))</span> |
---|
839 | <a name="ln-471"></a><span class="nl"> </span> <span class="n">DZ</span><span class="p">(</span><span class="n">LLH</span><span class="p">)</span><span class="o">=</span> |
---|
840 | <a name="ln-472"></a><span class="nl"> </span><span class="gs">*</span><span class="nb">DBLE</span><span class="p">(</span><span class="n">J</span><span class="p">)</span><span class="o">*</span><span class="n">WURZELB</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">-</span><span class="nb">DBLE</span><span class="p">(</span><span class="n">J</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="n">WURZELA</span><span class="o">*</span><span class="n">Z</span><span class="p">(</span><span class="n">LLP</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
841 | <a name="ln-473"></a><span class="nl"> </span><span class="gs"> </span> <span class="n">WURZELA</span><span class="o">=</span><span class="n">WURZELB</span> |
---|
842 | <a name="ln-474"></a><span class="nl"> </span> <span class="n">ENDIF</span> |
---|
843 | <a name="ln-475"></a><span class="nl"> </span> <span class="n">LLH</span><span class="o">=</span><span class="n">LLH</span><span class="o">+</span><span class="mi">1</span> |
---|
844 | <a name="ln-476"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">1</span> |
---|
845 | <a name="ln-477"></a><span class="nl">2 </span> <span class="k">CONTINUE</span> |
---|
846 | <a name="ln-478"></a><span class="nl"> </span> <span class="n">LLP</span><span class="o">=</span><span class="n">LLP</span><span class="o">+</span><span class="mi">1</span> |
---|
847 | <a name="ln-479"></a><span class="nl">1 CON</span><span class="gs">T</span><span class="n">INUE</span> |
---|
848 | <a name="ln-480"></a><span class="nl"> RETU</span><span class="gs">R</span><span class="n">N</span> |
---|
849 | <a name="ln-481"></a><span class="nl"> END </span><span class="gs">S</span><span class="n">UBROUTINE</span> <span class="n">DPLGND</span> |
---|
850 | <a name="ln-482"></a> |
---|
851 | <a name="ln-483"></a> |
---|
852 | <a name="ln-484"></a><span class="c">* Spectral Filter of Sardeshmukh and Hoskins (1984, MWR)</span> |
---|
853 | <a name="ln-485"></a><span class="c">* MM=Spectral truncation of field</span> |
---|
854 | <a name="ln-486"></a><span class="c">* MMAX= Spectral truncation of filter</span> |
---|
855 | <a name="ln-487"></a><span class="c">*</span> |
---|
856 | <a name="ln-488"></a><span class="nl"> </span> <span class="k">SUBROUTINE </span><span class="n">SPFILTER</span><span class="p">(</span><span class="n">FELDMN</span><span class="p">,</span><span class="n">MM</span><span class="p">,</span><span class="n">MMAX</span><span class="p">)</span> |
---|
857 | <a name="ln-489"></a> |
---|
858 | <a name="ln-490"></a><span class="nl"> </span> <span class="k">IMPLICIT NONE</span> |
---|
859 | <a name="ln-491"></a> |
---|
860 | <a name="ln-492"></a><span class="nl"> </span> <span class="kt">INTEGER </span><span class="n">MM</span><span class="p">,</span><span class="n">MMAX</span><span class="p">,</span><span class="n">I</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">K</span><span class="p">,</span><span class="n">L</span> |
---|
861 | <a name="ln-493"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">FELDMN</span><span class="p">(</span><span class="mi">0</span><span class="p">:(</span><span class="n">MM</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="p">(</span><span class="n">MM</span><span class="o">+</span><span class="mi">2</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> |
---|
862 | <a name="ln-494"></a><span class="nl"> </span> <span class="kt">REAL </span><span class="n">KMAX</span><span class="p">,</span><span class="n">SMAX</span><span class="p">,</span><span class="n">FAK</span> |
---|
863 | <a name="ln-495"></a> |
---|
864 | <a name="ln-496"></a><span class="nl"> </span> <span class="n">SMAX</span><span class="o">=</span><span class="mf">0.1</span> |
---|
865 | <a name="ln-497"></a><span class="nl"> </span> <span class="n">KMAX</span><span class="o">=-</span><span class="nb">ALOG</span><span class="p">(</span><span class="n">SMAX</span><span class="p">)</span> |
---|
866 | <a name="ln-498"></a><span class="nl"> </span> <span class="n">KMAX</span><span class="o">=</span><span class="n">KMAX</span><span class="o">/</span><span class="p">(</span><span class="nb">float</span><span class="p">(</span><span class="n">MMAX</span><span class="p">)</span><span class="o">*</span><span class="nb">float</span><span class="p">(</span><span class="n">MMAX</span><span class="o">+</span><span class="mi">1</span><span class="p">))</span><span class="o">**</span><span class="mi">2</span> |
---|
867 | <a name="ln-499"></a><span class="c">c WRITE(*,*)'alogsmax',alog(smax),'KMAX:',KMAX</span> |
---|
868 | <a name="ln-500"></a><span class="nl"> </span> <span class="n">l</span><span class="o">=</span><span class="mi">0</span> |
---|
869 | <a name="ln-501"></a><span class="nl"> </span> <span class="k">do </span><span class="n">i</span><span class="o">=</span><span class="mi">0</span><span class="p">,</span><span class="n">MM</span> |
---|
870 | <a name="ln-502"></a><span class="nl"> </span> <span class="k">do </span><span class="n">j</span><span class="o">=</span><span class="n">i</span><span class="p">,</span><span class="n">MM</span> |
---|
871 | <a name="ln-503"></a><span class="c">c write(*,*) i,j,feld(k),feld(k)*exp(-KMAX*(j*(j+1))**2)</span> |
---|
872 | <a name="ln-504"></a><span class="nl"> </span> <span class="k">if</span><span class="p">(</span><span class="n">j</span> <span class="p">.</span><span class="n">le</span><span class="p">.</span> <span class="n">MMAX</span><span class="p">)</span> <span class="k">then</span> |
---|
873 | <a name="ln-505"></a><span class="c">c fak=exp(-KMAX*(j*(j+1))**2)</span> |
---|
874 | <a name="ln-506"></a><span class="nl"> </span> <span class="n">fak</span><span class="o">=</span><span class="mf">1.0</span> |
---|
875 | <a name="ln-507"></a><span class="nl"> </span> <span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="p">)</span><span class="o">=</span><span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="p">)</span><span class="o">*</span><span class="n">fak</span> |
---|
876 | <a name="ln-508"></a><span class="nl"> </span> <span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="n">fak</span> |
---|
877 | <a name="ln-509"></a><span class="nl"> </span> <span class="k">else</span> |
---|
878 | <a name="ln-510"></a><span class="nl"> </span> <span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="p">)</span><span class="o">=</span><span class="mf">0.</span> |
---|
879 | <a name="ln-511"></a><span class="nl"> </span> <span class="n">feldmn</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">l</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span><span class="o">=</span><span class="mf">0.</span> |
---|
880 | <a name="ln-512"></a><span class="nl"> </span> <span class="n">endif</span> |
---|
881 | <a name="ln-513"></a><span class="nl"> </span> <span class="n">l</span><span class="o">=</span><span class="n">l</span><span class="o">+</span><span class="mi">1</span> |
---|
882 | <a name="ln-514"></a><span class="nl"> </span> <span class="n">enddo</span> |
---|
883 | <a name="ln-515"></a><span class="nl"> </span> <span class="n">enddo</span> |
---|
884 | <a name="ln-516"></a><span class="nl"> </span> <span class="k">END SUBROUTINE </span><span class="n">SPFILTER</span> |
---|
885 | <a name="ln-517"></a> |
---|
886 | <a name="ln-518"></a><span class="nl"> </span> <span class="k">END MODULE </span><span class="n">PHTOGR</span> |
---|
887 | <a name="ln-519"></a> |
---|
888 | </pre></div> |
---|
889 | |
---|
890 | </section> |
---|
891 | </div> |
---|
892 | </div> |
---|
893 | |
---|
894 | <section class="visible-xs visible-sm hidden-md"> |
---|
895 | <hr> |
---|
896 | |
---|
897 | |
---|
898 | <div class="panel panel-default"> |
---|
899 | <div class="panel-heading text-left"><h3 class="panel-title"><a data-toggle="collapse" href="#allfiles-1">All Source Files</a></h3></div> |
---|
900 | <div id="allfiles-1" class="panel-collapse collapse"> |
---|
901 | <div class="list-group"> |
---|
902 | |
---|
903 | <a class="list-group-item" href="../sourcefile/grphreal.f.html">grphreal.f</a> |
---|
904 | |
---|
905 | <a class="list-group-item" href="../sourcefile/phgrreal.f.html">phgrreal.f</a> |
---|
906 | |
---|
907 | <a class="list-group-item" href="../sourcefile/preconvert.f90.html">preconvert.f90</a> |
---|
908 | |
---|
909 | <a class="list-group-item" href="../sourcefile/rwgrib2.f90.html">rwGRIB2.f90</a> |
---|
910 | |
---|
911 | </div> |
---|
912 | </div> |
---|
913 | </div> |
---|
914 | |
---|
915 | |
---|
916 | </section> |
---|
917 | |
---|
918 | <hr> |
---|
919 | </div> <!-- /container --> |
---|
920 | <footer> |
---|
921 | <div class="container"> |
---|
922 | <div class="row"> |
---|
923 | <div class="col-xs-6 col-md-4"><p>© 2019 <a rel="license" href="http://creativecommons.org/licenses/by/4.0/"><img alt="Creative Commons License" style="border-width:0" src="../80x15.png" /></a></p></div> |
---|
924 | <div class="col-xs-6 col-md-4 col-md-push-4"> |
---|
925 | <p class="text-right"> |
---|
926 | Documentation generated by |
---|
927 | <a href="https://github.com/cmacmackin/ford">FORD</a> |
---|
928 | |
---|
929 | </p> |
---|
930 | </div> |
---|
931 | <div class="col-xs-12 col-md-4 col-md-pull-4"><p class="text-center"> Flex_extract: Calculation of etadot was developed by Leopold Haimberger</p></div> |
---|
932 | </div> |
---|
933 | <br> |
---|
934 | </div> <!-- /container --> |
---|
935 | </footer> |
---|
936 | |
---|
937 | <!-- Bootstrap core JavaScript |
---|
938 | ================================================== --> |
---|
939 | <!-- Placed at the end of the document so the pages load faster --> |
---|
940 | <!-- |
---|
941 | <script src="https://ajax.googleapis.com/ajax/libs/jquery/1.11.1/jquery.min.js"></script> |
---|
942 | --> |
---|
943 | <script src="../js/bootstrap.min.js"></script> |
---|
944 | <!-- IE10 viewport hack for Surface/desktop Windows 8 bug --> |
---|
945 | <script src="../js/ie10-viewport-bug-workaround.js"></script> |
---|
946 | |
---|
947 | <!-- MathJax JavaScript |
---|
948 | ================================================== --> |
---|
949 | <!-- Placed at the end of the document so the pages load faster --> |
---|
950 | <script type="text/x-mathjax-config"> |
---|
951 | MathJax.Hub.Config({ |
---|
952 | TeX: { extensions: ['AMSmath.js','AMSsymbols.js','noErrors.js','noUndefined.js'], equationNumbers: { autoNumber: 'AMS' } }, |
---|
953 | jax: ['input/TeX','input/MathML','output/HTML-CSS'], |
---|
954 | extensions: ['tex2jax.js','mml2jax.js','MathMenu.js','MathZoom.js'], |
---|
955 | 'HTML-CSS': { |
---|
956 | styles: { '.MathJax_Display, .MathJax .mo, .MathJax .mi, .MathJax .mn': {color: '#000000 ! important'} } |
---|
957 | } |
---|
958 | }); |
---|
959 | </script> |
---|
960 | <!-- <script src="https://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML"></script> --> |
---|
961 | </body> |
---|
962 | </html> |
---|