legsym_8f90_source.html 49 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700
  1. <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
  2. <html xmlns="http://www.w3.org/1999/xhtml">
  3. <head>
  4. <meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
  5. <title>PUMA: /Users/home/WC/puma/src/legsym.f90 Source File</title>
  6. <link href="tabs.css" rel="stylesheet" type="text/css"/>
  7. <link href="doxygen.css" rel="stylesheet" type="text/css" />
  8. <link href="navtree.css" rel="stylesheet" type="text/css"/>
  9. <script type="text/javascript" src="jquery.js"></script>
  10. <script type="text/javascript" src="resize.js"></script>
  11. <script type="text/javascript" src="navtree.js"></script>
  12. <script type="text/javascript">
  13. $(document).ready(initResizable);
  14. </script>
  15. <link href="search/search.css" rel="stylesheet" type="text/css"/>
  16. <script type="text/javascript" src="search/search.js"></script>
  17. <script type="text/javascript">
  18. $(document).ready(function() { searchBox.OnSelectItem(0); });
  19. </script>
  20. </head>
  21. <body>
  22. <div id="top"><!-- do not remove this div! -->
  23. <div id="titlearea">
  24. <table cellspacing="0" cellpadding="0">
  25. <tbody>
  26. <tr style="height: 56px;">
  27. <td id="projectlogo"><img alt="Logo" src="puma103.jpg"/></td>
  28. <td style="padding-left: 0.5em;">
  29. <div id="projectname">PUMA
  30. &#160;<span id="projectnumber">219</span>
  31. </div>
  32. <div id="projectbrief">Portable University Model of the Atmosphere</div>
  33. </td>
  34. </tr>
  35. </tbody>
  36. </table>
  37. </div>
  38. <!-- Generated by Doxygen 1.7.5.1 -->
  39. <script type="text/javascript">
  40. var searchBox = new SearchBox("searchBox", "search",false,'Search');
  41. </script>
  42. <div id="navrow1" class="tabs">
  43. <ul class="tablist">
  44. <li><a href="index.html"><span>Main&#160;Page</span></a></li>
  45. <li><a href="annotated.html"><span>Data&#160;Types&#160;List</span></a></li>
  46. <li class="current"><a href="files.html"><span>Files</span></a></li>
  47. <li>
  48. <div id="MSearchBox" class="MSearchBoxInactive">
  49. <span class="left">
  50. <img id="MSearchSelect" src="search/mag_sel.png"
  51. onmouseover="return searchBox.OnSearchSelectShow()"
  52. onmouseout="return searchBox.OnSearchSelectHide()"
  53. alt=""/>
  54. <input type="text" id="MSearchField" value="Search" accesskey="S"
  55. onfocus="searchBox.OnSearchFieldFocus(true)"
  56. onblur="searchBox.OnSearchFieldFocus(false)"
  57. onkeyup="searchBox.OnSearchFieldChange(event)"/>
  58. </span><span class="right">
  59. <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
  60. </span>
  61. </div>
  62. </li>
  63. </ul>
  64. </div>
  65. <div id="navrow2" class="tabs2">
  66. <ul class="tablist">
  67. <li><a href="files.html"><span>File&#160;List</span></a></li>
  68. <li><a href="globals.html"><span>File&#160;Members</span></a></li>
  69. </ul>
  70. </div>
  71. </div>
  72. <div id="side-nav" class="ui-resizable side-nav-resizable">
  73. <div id="nav-tree">
  74. <div id="nav-tree-contents">
  75. </div>
  76. </div>
  77. <div id="splitbar" style="-moz-user-select:none;"
  78. class="ui-resizable-handle">
  79. </div>
  80. </div>
  81. <script type="text/javascript">
  82. initNavTree('legsym_8f90.html','');
  83. </script>
  84. <div id="doc-content">
  85. <div class="header">
  86. <div class="headertitle">
  87. <div class="title">/Users/home/WC/puma/src/legsym.f90</div> </div>
  88. </div>
  89. <div class="contents">
  90. <a href="legsym_8f90.html">Go to the documentation of this file.</a><div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 <span class="comment">! ************************************************************************</span>
  91. <a name="l00002"></a>00002 <span class="comment">! * module legsym - direct and indirect Legendre transformation routines *</span>
  92. <a name="l00003"></a>00003 <span class="comment">! * using symmetric and antisymmetric fourier coefficients *</span>
  93. <a name="l00004"></a>00004 <span class="comment">! * E. Kirk 08-Sep-2010 tested for T21 - T682 resolutions 32 &amp; 64 bit *</span>
  94. <a name="l00005"></a>00005 <span class="comment">! ************************************************************************</span>
  95. <a name="l00006"></a>00006
  96. <a name="l00007"></a><a class="code" href="classlegsym.html">00007</a> <span class="keyword">module</span> <a class="code" href="classlegsym.html">legsym</a>
  97. <a name="l00008"></a>00008 <span class="comment">! ************************</span>
  98. <a name="l00009"></a>00009 <span class="comment">! * Legendre Polynomials *</span>
  99. <a name="l00010"></a>00010 <span class="comment">! ************************</span>
  100. <a name="l00011"></a>00011
  101. <a name="l00012"></a><a class="code" href="classlegsym.html#a9f44e2174fc92267b46632e481acce37">00012</a> <span class="keywordtype">integer</span> :: ntru
  102. <a name="l00013"></a><a class="code" href="classlegsym.html#a3de2526b7367430bf2925e4d9e293b8d">00013</a> <span class="keywordtype">integer</span> :: ntp1
  103. <a name="l00014"></a><a class="code" href="classlegsym.html#a1cc07c984867c3dc3b122ac8fce446e2">00014</a> <span class="keywordtype">integer</span> :: ncsp
  104. <a name="l00015"></a><a class="code" href="classlegsym.html#a8d1597f9974455b1abb4b57e34a77625">00015</a> <span class="keywordtype">integer</span> :: nesp
  105. <a name="l00016"></a><a class="code" href="classlegsym.html#adc6b16f2626e9dcb618d2d950e2db17a">00016</a> <span class="keywordtype">integer</span> :: nlon
  106. <a name="l00017"></a><a class="code" href="classlegsym.html#a296b840ee1daf32583f78f3b8136acb0">00017</a> <span class="keywordtype">integer</span> :: nlpp
  107. <a name="l00018"></a><a class="code" href="classlegsym.html#a2ca038107ed480d5443db9a2f0322ce1">00018</a> <span class="keywordtype">integer</span> :: nhpp
  108. <a name="l00019"></a><a class="code" href="classlegsym.html#a48c28c34fad8572472b4994bf99da184">00019</a> <span class="keywordtype">integer</span> :: nlat
  109. <a name="l00020"></a><a class="code" href="classlegsym.html#a8b4c67467362c9d324d5ad7bf5b6e78e">00020</a> <span class="keywordtype">integer</span> :: nlev
  110. <a name="l00021"></a><a class="code" href="classlegsym.html#adb777017fc0833c5be3bbe40fb789d9f">00021</a> <span class="keywordtype">real</span> :: plavor
  111. <a name="l00022"></a>00022
  112. <a name="l00023"></a><a class="code" href="classlegsym.html#a35897a5a32d0ac4f7193e2c6df059430">00023</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qi(:,:) <span class="comment">! P(m,n) = Associated Legendre Polynomials</span>
  113. <a name="l00024"></a><a class="code" href="classlegsym.html#a3fd3514e19855e9c94debeb13d1c0a3c">00024</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qj(:,:) <span class="comment">! Q(m,n) = Used for d/d(mu)</span>
  114. <a name="l00025"></a><a class="code" href="classlegsym.html#a8da266cccc13b41ffb15d593f3cfbda1">00025</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qc(:,:) <span class="comment">! P(m,n) * gwd used in fc2sp</span>
  115. <a name="l00026"></a><a class="code" href="classlegsym.html#a27a404e273f903a7a9fa790dd2440605">00026</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qe(:,:) <span class="comment">! Q(mn,) * gwd / cos2 used in mktend</span>
  116. <a name="l00027"></a><a class="code" href="classlegsym.html#ae957dbe2cb8c53bef3e92d64bbea2785">00027</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qq(:,:) <span class="comment">! P(m,n) * gwd / cos2 * n * (n+1) / 2 &quot;</span>
  117. <a name="l00028"></a><a class="code" href="classlegsym.html#a9e9b5d567f315f6cf660a0b44ffb4889">00028</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qu(:,:) <span class="comment">! P(m,n) / (n*(n+1)) * m used in dv2uv</span>
  118. <a name="l00029"></a><a class="code" href="classlegsym.html#a9a8f54ece20fd3aacb5d1f1570450ee0">00029</a> <span class="keywordtype">real</span> , <span class="keywordtype">allocatable</span> :: qv(:,:) <span class="comment">! Q(m,n) / (n*(n+1)) used in dv2uv</span>
  119. <a name="l00030"></a><a class="code" href="classlegsym.html#aff59025facb5c5b3d3717ca32328f587">00030</a> <span class="keywordtype">complex</span>, <span class="keywordtype">allocatable</span> :: qx(:,:) <span class="comment">! P(m,n) * gwd / cos2 * m used in mktend</span>
  120. <a name="l00031"></a>00031 <span class="keyword">end module legsym</span>
  121. <a name="l00032"></a>00032
  122. <a name="l00033"></a>00033 <span class="comment">! =================</span>
  123. <a name="l00034"></a>00034 <span class="comment">! SUBROUTINE LEGINI</span>
  124. <a name="l00035"></a>00035 <span class="comment">! =================</span>
  125. <a name="l00036"></a>00036
  126. <a name="l00037"></a><a class="code" href="legsym_8f90.html#a86bc436e65d6c4ddde72bb3cce7dc8c8">00037</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a86bc436e65d6c4ddde72bb3cce7dc8c8">legini</a>(klat,klpp,kesp,klev,vorpla,sid,gwd)
  127. <a name="l00038"></a>00038 use <span class="keywordflow">legsym</span>
  128. <a name="l00039"></a>00039 <span class="keyword">implicit none</span>
  129. <a name="l00040"></a>00040
  130. <a name="l00041"></a>00041 <span class="keywordtype">integer</span> :: klat
  131. <a name="l00042"></a>00042 <span class="keywordtype">integer</span> :: klpp
  132. <a name="l00043"></a>00043 <span class="keywordtype">integer</span> :: kesp
  133. <a name="l00044"></a>00044 <span class="keywordtype">integer</span> :: klev
  134. <a name="l00045"></a>00045
  135. <a name="l00046"></a>00046 <span class="keywordtype">real</span> :: vorpla <span class="comment">! planetary vorticity</span>
  136. <a name="l00047"></a>00047 <span class="keywordtype">real (kind=8)</span> :: sid(*) <span class="comment">! sin(phi)</span>
  137. <a name="l00048"></a>00048 <span class="keywordtype">real (kind=8)</span> :: gwd(*) <span class="comment">! Gaussian weight (phi)</span>
  138. <a name="l00049"></a>00049
  139. <a name="l00050"></a>00050 <span class="keywordtype">integer</span> :: jlat <span class="comment">! Latitude</span>
  140. <a name="l00051"></a>00051 <span class="keywordtype">integer</span> :: lm
  141. <a name="l00052"></a>00052 <span class="keywordtype">integer</span> :: m
  142. <a name="l00053"></a>00053 <span class="keywordtype">integer</span> :: n
  143. <a name="l00054"></a>00054
  144. <a name="l00055"></a>00055 <span class="keywordtype">real (kind=8)</span> :: amsq
  145. <a name="l00056"></a>00056 <span class="keywordtype">real (kind=8)</span> :: z1
  146. <a name="l00057"></a>00057 <span class="keywordtype">real (kind=8)</span> :: z2
  147. <a name="l00058"></a>00058 <span class="keywordtype">real (kind=8)</span> :: z3
  148. <a name="l00059"></a>00059 <span class="keywordtype">real (kind=8)</span> :: f1m
  149. <a name="l00060"></a>00060 <span class="keywordtype">real (kind=8)</span> :: f2m
  150. <a name="l00061"></a>00061 <span class="keywordtype">real (kind=8)</span> :: znn1
  151. <a name="l00062"></a>00062 <span class="keywordtype">real (kind=8)</span> :: zsin <span class="comment">! sin</span>
  152. <a name="l00063"></a>00063 <span class="keywordtype">real (kind=8)</span> :: zcsq <span class="comment">! cos2</span>
  153. <a name="l00064"></a>00064 <span class="keywordtype">real (kind=8)</span> :: zcos <span class="comment">! cos</span>
  154. <a name="l00065"></a>00065 <span class="keywordtype">real (kind=8)</span> :: zgwd <span class="comment">! gw</span>
  155. <a name="l00066"></a>00066 <span class="keywordtype">real (kind=8)</span> :: zgwdcsq <span class="comment">! gw / cos2</span>
  156. <a name="l00067"></a>00067 <span class="keywordtype">real (kind=8)</span> :: zpli(kesp)
  157. <a name="l00068"></a>00068 <span class="keywordtype">real (kind=8)</span> :: zpld(kesp)
  158. <a name="l00069"></a>00069
  159. <a name="l00070"></a>00070 nlat = klat
  160. <a name="l00071"></a>00071 nlpp = klpp
  161. <a name="l00072"></a>00072 nhpp = klpp / 2
  162. <a name="l00073"></a>00073 nlon = klat + klat
  163. <a name="l00074"></a>00074 ntru = (nlon - 1) /3
  164. <a name="l00075"></a>00075 ntp1 = ntru + 1
  165. <a name="l00076"></a>00076 ncsp = ((ntru + 1) * (ntru + 2)) / 2
  166. <a name="l00077"></a>00077 nesp = kesp
  167. <a name="l00078"></a>00078 nlev = klev
  168. <a name="l00079"></a>00079
  169. <a name="l00080"></a>00080 plavor = vorpla
  170. <a name="l00081"></a>00081
  171. <a name="l00082"></a>00082 <span class="keyword">allocate</span>(qi(ncsp,nhpp))
  172. <a name="l00083"></a>00083 <span class="keyword">allocate</span>(qj(ncsp,nhpp))
  173. <a name="l00084"></a>00084 <span class="keyword">allocate</span>(qc(ncsp,nhpp))
  174. <a name="l00085"></a>00085 <span class="keyword">allocate</span>(qe(ncsp,nhpp))
  175. <a name="l00086"></a>00086 <span class="keyword">allocate</span>(qx(ncsp,nhpp))
  176. <a name="l00087"></a>00087 <span class="keyword">allocate</span>(qq(ncsp,nhpp))
  177. <a name="l00088"></a>00088 <span class="keyword">allocate</span>(qu(ncsp,nhpp))
  178. <a name="l00089"></a>00089 <span class="keyword">allocate</span>(qv(ncsp,nhpp))
  179. <a name="l00090"></a>00090
  180. <a name="l00091"></a>00091 <span class="keyword">do</span> jlat = 1 , nhpp
  181. <a name="l00092"></a>00092
  182. <a name="l00093"></a>00093 <span class="comment">! set p(0,0) and p(0,1)</span>
  183. <a name="l00094"></a>00094
  184. <a name="l00095"></a>00095 zgwd = gwd(jlat) <span class="comment">! gaussian weight - from inigau</span>
  185. <a name="l00096"></a>00096 zsin = sid(jlat) <span class="comment">! sin(phi) - from inigau</span>
  186. <a name="l00097"></a>00097 zcsq = 1.0_8 - zsin * zsin <span class="comment">! cos(phi) squared</span>
  187. <a name="l00098"></a>00098 zgwdcsq = zgwd / zcsq <span class="comment">! weight / cos squared</span>
  188. <a name="l00099"></a>00099 zcos = sqrt(zcsq) <span class="comment">! cos(phi)</span>
  189. <a name="l00100"></a>00100 f1m = sqrt(1.5_8)
  190. <a name="l00101"></a>00101 zpli(1) = sqrt(0.5_8)
  191. <a name="l00102"></a>00102 zpli(2) = f1m * zsin
  192. <a name="l00103"></a>00103 zpld(1) = 0.0
  193. <a name="l00104"></a>00104 lm = 2
  194. <a name="l00105"></a>00105
  195. <a name="l00106"></a>00106 <span class="comment">! loop over wavenumbers</span>
  196. <a name="l00107"></a>00107
  197. <a name="l00108"></a>00108 <span class="keyword">do</span> m = 0 , ntru
  198. <a name="l00109"></a>00109 <span class="keyword">if</span> (m &gt; 0) <span class="keyword">then</span>
  199. <a name="l00110"></a>00110 lm = lm + 1
  200. <a name="l00111"></a>00111 f2m = -f1m * sqrt(zcsq / (m+m))
  201. <a name="l00112"></a>00112 f1m = f2m * sqrt(m+m + 3.0_8)
  202. <a name="l00113"></a>00113 zpli(lm) = f2m
  203. <a name="l00114"></a>00114 <span class="keyword">if</span> (lm &lt; ncsp) <span class="keyword">then</span>
  204. <a name="l00115"></a>00115 lm = lm + 1
  205. <a name="l00116"></a>00116 zpli(lm ) = f1m * zsin
  206. <a name="l00117"></a>00117 zpld(lm-1) = -m * f2m * zsin
  207. <a name="l00118"></a>00118 <span class="keyword">endif</span> <span class="comment">! (lm &lt; ncsp)</span>
  208. <a name="l00119"></a>00119 <span class="keyword">endif</span> <span class="comment">! (m &gt; 0)</span>
  209. <a name="l00120"></a>00120
  210. <a name="l00121"></a>00121 amsq = m * m
  211. <a name="l00122"></a>00122
  212. <a name="l00123"></a>00123 <span class="keyword">do</span> n = m+2 , ntru
  213. <a name="l00124"></a>00124 lm = lm + 1
  214. <a name="l00125"></a>00125 z1 = sqrt(((n-1)*(n-1) - amsq) / (4*(n-1)*(n-1)-1))
  215. <a name="l00126"></a>00126 z2 = zsin * zpli(lm-1) - z1 * zpli(lm-2)
  216. <a name="l00127"></a>00127 zpli(lm ) = z2 * sqrt((4*n*n-1) / (n*n-amsq))
  217. <a name="l00128"></a>00128 zpld(lm-1) = (1-n) * z2 + n * z1 * zpli(lm-2)
  218. <a name="l00129"></a>00129 <span class="keyword">enddo</span> <span class="comment">! n</span>
  219. <a name="l00130"></a>00130
  220. <a name="l00131"></a>00131 <span class="keyword">if</span> (lm &lt; ncsp) <span class="keyword">then</span> <span class="comment">! mode (m,ntru)</span>
  221. <a name="l00132"></a>00132 z3 = sqrt((ntru*ntru-amsq) / (4*ntru*ntru-1))
  222. <a name="l00133"></a>00133 zpld(lm)=-ntru*zsin*zpli(lm) + (ntru+ntru+1)*zpli(lm-1)*z3
  223. <a name="l00134"></a>00134 <span class="keyword">else</span> <span class="comment">! mode (ntru,ntru)</span>
  224. <a name="l00135"></a>00135 zpld(lm)=-ntru*zsin*zpli(lm)
  225. <a name="l00136"></a>00136 <span class="keyword">endif</span>
  226. <a name="l00137"></a>00137 <span class="keyword">enddo</span> <span class="comment">! m</span>
  227. <a name="l00138"></a>00138
  228. <a name="l00139"></a>00139 lm = 0
  229. <a name="l00140"></a>00140 <span class="keyword">do</span> m = 0 , ntru
  230. <a name="l00141"></a>00141 <span class="keyword">do</span> n = m , ntru
  231. <a name="l00142"></a>00142 lm = lm + 1
  232. <a name="l00143"></a>00143 znn1 = 0.0
  233. <a name="l00144"></a>00144 <span class="keyword">if</span> (n &gt; 0) znn1 = 1.0_8 / (n*(n+1))
  234. <a name="l00145"></a>00145 qi(lm,jlat) = zpli(lm)
  235. <a name="l00146"></a>00146 qj(lm,jlat) = zpld(lm)
  236. <a name="l00147"></a>00147 qc(lm,jlat) = zpli(lm) * zgwd
  237. <a name="l00148"></a>00148 qu(lm,jlat) = zpli(lm) * znn1 * m
  238. <a name="l00149"></a>00149 qv(lm,jlat) = zpld(lm) * znn1
  239. <a name="l00150"></a>00150 qe(lm,jlat) = zpld(lm) * zgwdcsq
  240. <a name="l00151"></a>00151 qq(lm,jlat) = zpli(lm) * zgwdcsq * n * (n+1) * 0.5_8
  241. <a name="l00152"></a>00152 qx(lm,jlat) = zpli(lm) * zgwdcsq * m * (0.0,1.0)
  242. <a name="l00153"></a>00153 <span class="keyword">enddo</span> <span class="comment">! n</span>
  243. <a name="l00154"></a>00154 <span class="keyword">enddo</span> <span class="comment">! m</span>
  244. <a name="l00155"></a>00155 <span class="keyword">enddo</span> <span class="comment">! jlat</span>
  245. <a name="l00156"></a>00156 return
  246. <a name="l00157"></a>00157 <span class="keyword">end</span>
  247. <a name="l00158"></a>00158
  248. <a name="l00159"></a>00159
  249. <a name="l00160"></a>00160 <span class="comment">! ================</span>
  250. <a name="l00161"></a>00161 <span class="comment">! SUBROUTINE FC2SP</span>
  251. <a name="l00162"></a>00162 <span class="comment">! ================</span>
  252. <a name="l00163"></a>00163
  253. <a name="l00164"></a><a class="code" href="legsym_8f90.html#a04d46a94caf6c743547ac25cfa3058d4">00164</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a04d46a94caf6c743547ac25cfa3058d4">fc2sp</a>(fc,sp)
  254. <a name="l00165"></a>00165 use <span class="keywordflow">legsym</span>
  255. <a name="l00166"></a>00166 <span class="keyword">implicit none</span>
  256. <a name="l00167"></a>00167 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in )</span> :: fc(nlon,nhpp)
  257. <a name="l00168"></a>00168 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(out)</span> :: sp(nesp/2)
  258. <a name="l00169"></a>00169
  259. <a name="l00170"></a>00170 <span class="keywordtype">integer</span> :: l <span class="comment">! Index for latitude</span>
  260. <a name="l00171"></a>00171 <span class="keywordtype">integer</span> :: m <span class="comment">! Index for zonal wavenumber</span>
  261. <a name="l00172"></a>00172 <span class="keywordtype">integer</span> :: w <span class="comment">! Index for spherical harmonic</span>
  262. <a name="l00173"></a>00173 <span class="keywordtype">integer</span> :: e <span class="comment">! Index for last wavenumber</span>
  263. <a name="l00174"></a>00174
  264. <a name="l00175"></a>00175 sp(:) = (0.0,0.0)
  265. <a name="l00176"></a>00176
  266. <a name="l00177"></a>00177 <span class="keyword">do</span> l = 1 , nhpp
  267. <a name="l00178"></a>00178 w = 1
  268. <a name="l00179"></a>00179 <span class="keyword">do</span> m = 1 , ntp1
  269. <a name="l00180"></a>00180 e = w + ntp1 - m
  270. <a name="l00181"></a>00181 sp(w :e:2) = sp(w :e:2) + qc(w :e:2,l) * (fc(m,l) + fc(m+nlat,l))
  271. <a name="l00182"></a>00182 sp(w+1:e:2) = sp(w+1:e:2) + qc(w+1:e:2,l) * (fc(m,l) - fc(m+nlat,l))
  272. <a name="l00183"></a>00183 w = e + 1
  273. <a name="l00184"></a>00184 <span class="keyword">enddo</span> <span class="comment">! m</span>
  274. <a name="l00185"></a>00185 <span class="keyword">enddo</span> <span class="comment">! l</span>
  275. <a name="l00186"></a>00186 return
  276. <a name="l00187"></a>00187 <span class="keyword">end</span>
  277. <a name="l00188"></a>00188
  278. <a name="l00189"></a>00189
  279. <a name="l00190"></a>00190 <span class="comment">! ================</span>
  280. <a name="l00191"></a>00191 <span class="comment">! SUBROUTINE SP2FC</span>
  281. <a name="l00192"></a>00192 <span class="comment">! ================</span>
  282. <a name="l00193"></a>00193
  283. <a name="l00194"></a><a class="code" href="legsym_8f90.html#aec404fc15930c6e4584c088a399ea099">00194</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#aec404fc15930c6e4584c088a399ea099">sp2fc</a>(sp,fc) ! Spectral to Fourier
  284. <a name="l00195"></a>00195 use <span class="keywordflow">legsym</span>
  285. <a name="l00196"></a>00196 <span class="keyword">implicit none</span>
  286. <a name="l00197"></a>00197
  287. <a name="l00198"></a>00198 <span class="keywordtype">complex</span> :: sp(ncsp) <span class="comment">! Coefficients of spherical harmonics</span>
  288. <a name="l00199"></a>00199 <span class="keywordtype">complex</span> :: fc(nlon,nhpp) <span class="comment">! Fourier coefficients</span>
  289. <a name="l00200"></a>00200
  290. <a name="l00201"></a>00201 <span class="keywordtype">integer</span> :: l <span class="comment">! Loop index for latitude</span>
  291. <a name="l00202"></a>00202 <span class="keywordtype">integer</span> :: m <span class="comment">! Loop index for zonal wavenumber m</span>
  292. <a name="l00203"></a>00203 <span class="keywordtype">integer</span> :: w <span class="comment">! Index for spectral mode</span>
  293. <a name="l00204"></a>00204 <span class="keywordtype">integer</span> :: e <span class="comment">! Index for last wavenumber</span>
  294. <a name="l00205"></a>00205
  295. <a name="l00206"></a>00206 <span class="keywordtype">complex</span> :: fs,fa
  296. <a name="l00207"></a>00207
  297. <a name="l00208"></a>00208 fc(:,:) = (0.0,0.0)
  298. <a name="l00209"></a>00209
  299. <a name="l00210"></a>00210 <span class="keyword">do</span> l = 1 , nhpp
  300. <a name="l00211"></a>00211 w = 1
  301. <a name="l00212"></a>00212 <span class="keyword">do</span> m = 1 ,ntp1
  302. <a name="l00213"></a>00213 e = w + ntp1 - m
  303. <a name="l00214"></a>00214 fs = dot_product(qi(w :e:2,l),sp(w :e:2))
  304. <a name="l00215"></a>00215 fa = dot_product(qi(w+1:e:2,l),sp(w+1:e:2))
  305. <a name="l00216"></a>00216 fc(m ,l) = fs + fa
  306. <a name="l00217"></a>00217 fc(m+nlat,l) = fs - fa
  307. <a name="l00218"></a>00218 w = e + 1
  308. <a name="l00219"></a>00219 <span class="keyword">enddo</span> <span class="comment">! m</span>
  309. <a name="l00220"></a>00220 <span class="keyword">enddo</span> <span class="comment">! l</span>
  310. <a name="l00221"></a>00221 return
  311. <a name="l00222"></a>00222 <span class="keyword">end</span>
  312. <a name="l00223"></a>00223
  313. <a name="l00224"></a>00224
  314. <a name="l00225"></a>00225 <span class="comment">! ===================</span>
  315. <a name="l00226"></a>00226 <span class="comment">! SUBROUTINE SP2FCDMU</span>
  316. <a name="l00227"></a>00227 <span class="comment">! ===================</span>
  317. <a name="l00228"></a>00228
  318. <a name="l00229"></a><a class="code" href="legsym_8f90.html#ac25a3c42ee19118b299203d2747cb59e">00229</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#ac25a3c42ee19118b299203d2747cb59e">sp2fcdmu</a>(sp,fc) ! Spectral to Fourier d/dmu
  319. <a name="l00230"></a>00230 use <span class="keywordflow">legsym</span>
  320. <a name="l00231"></a>00231 <span class="keyword">implicit none</span>
  321. <a name="l00232"></a>00232
  322. <a name="l00233"></a>00233 <span class="keywordtype">complex</span> :: sp(ncsp) <span class="comment">! Coefficients of spherical harmonics</span>
  323. <a name="l00234"></a>00234 <span class="keywordtype">complex</span> :: fc(nlon,nhpp) <span class="comment">! Fourier coefficients</span>
  324. <a name="l00235"></a>00235
  325. <a name="l00236"></a>00236 <span class="keywordtype">integer</span> :: l <span class="comment">! Loop index for latitude</span>
  326. <a name="l00237"></a>00237 <span class="keywordtype">integer</span> :: m <span class="comment">! Loop index for zonal wavenumber m</span>
  327. <a name="l00238"></a>00238 <span class="keywordtype">integer</span> :: w <span class="comment">! Index for spectral mode</span>
  328. <a name="l00239"></a>00239 <span class="keywordtype">integer</span> :: e <span class="comment">! Index for last wavenumber</span>
  329. <a name="l00240"></a>00240
  330. <a name="l00241"></a>00241 <span class="keywordtype">complex</span> :: fs,fa
  331. <a name="l00242"></a>00242
  332. <a name="l00243"></a>00243 fc(:,:) = (0.0,0.0)
  333. <a name="l00244"></a>00244
  334. <a name="l00245"></a>00245 <span class="keyword">do</span> l = 1 , nhpp
  335. <a name="l00246"></a>00246 w = 1
  336. <a name="l00247"></a>00247 <span class="keyword">do</span> m = 1 , ntp1
  337. <a name="l00248"></a>00248 e = w + ntp1 - m
  338. <a name="l00249"></a>00249 fs = dot_product(qj(w :e:2,l),sp(w :e:2))
  339. <a name="l00250"></a>00250 fa = dot_product(qj(w+1:e:2,l),sp(w+1:e:2))
  340. <a name="l00251"></a>00251 fc(m ,l) = fa + fs
  341. <a name="l00252"></a>00252 fc(m+nlat,l) = fa - fs
  342. <a name="l00253"></a>00253 w = e + 1
  343. <a name="l00254"></a>00254 <span class="keyword">enddo</span> <span class="comment">! m</span>
  344. <a name="l00255"></a>00255 <span class="keyword">enddo</span> <span class="comment">! l</span>
  345. <a name="l00256"></a>00256 return
  346. <a name="l00257"></a>00257 <span class="keyword">end</span>
  347. <a name="l00258"></a>00258
  348. <a name="l00259"></a>00259
  349. <a name="l00260"></a>00260 <span class="comment">! ================</span>
  350. <a name="l00261"></a>00261 <span class="comment">! SUBROUTINE DV2UV</span>
  351. <a name="l00262"></a>00262 <span class="comment">! ================</span>
  352. <a name="l00263"></a>00263
  353. <a name="l00264"></a>00264 <span class="comment">! This is an alternative subroutine for computing U and V from Div and Vor.</span>
  354. <a name="l00265"></a>00265 <span class="comment">! It looks much prettier than the regular one, but unfortunately it is slower</span>
  355. <a name="l00266"></a>00266 <span class="comment">! if compiling with &quot;gfortran&quot;. </span>
  356. <a name="l00267"></a>00267 <span class="comment">! I leave it (unused) in this module for educational purposes.</span>
  357. <a name="l00268"></a>00268
  358. <a name="l00269"></a><a class="code" href="legsym_8f90.html#a581747beee32c3386dfc0b9b1cfa4e0f">00269</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a581747beee32c3386dfc0b9b1cfa4e0f">dv2uv_alt</a>(pd,pz,pu,pv)
  359. <a name="l00270"></a>00270 use <span class="keywordflow">legsym</span>
  360. <a name="l00271"></a>00271 <span class="keyword">implicit none</span>
  361. <a name="l00272"></a>00272 <span class="keywordtype">complex</span>, <span class="keywordtype">parameter</span> :: i = (0.0,1.0)
  362. <a name="l00273"></a>00273 <span class="keywordtype">complex</span> :: pd(nesp/2) <span class="comment">! Spherical harmonics of divergence</span>
  363. <a name="l00274"></a>00274 <span class="keywordtype">complex</span> :: pz(nesp/2) <span class="comment">! Spherical harmonics of vorticity</span>
  364. <a name="l00275"></a>00275 <span class="keywordtype">complex</span> :: pu(nlon,nhpp) <span class="comment">! Fourier coefficients of u</span>
  365. <a name="l00276"></a>00276 <span class="keywordtype">complex</span> :: pv(nlon,nhpp) <span class="comment">! Fourier coefficients of v</span>
  366. <a name="l00277"></a>00277 <span class="keywordtype">complex</span> :: zsave,uds,vds,uzs,vzs,uda,vda,uza,vza
  367. <a name="l00278"></a>00278
  368. <a name="l00279"></a>00279 <span class="keywordtype">integer</span> :: l <span class="comment">! Loop index for latitude</span>
  369. <a name="l00280"></a>00280 <span class="keywordtype">integer</span> :: m <span class="comment">! Loop index for zonal wavenumber m</span>
  370. <a name="l00281"></a>00281 <span class="keywordtype">integer</span> :: w <span class="comment">! Loop index for spectral mode</span>
  371. <a name="l00282"></a>00282 <span class="keywordtype">integer</span> :: e <span class="comment">! End index</span>
  372. <a name="l00283"></a>00283
  373. <a name="l00284"></a>00284 pu(:,:) = (0.0,0.0)
  374. <a name="l00285"></a>00285 pv(:,:) = (0.0,0.0)
  375. <a name="l00286"></a>00286
  376. <a name="l00287"></a>00287 zsave = pz(2) <span class="comment">! Save mode(0,1) of vorticity</span>
  377. <a name="l00288"></a>00288 pz(2) = zsave - cmplx(plavor,0.0) <span class="comment">! Convert pz from absolute to relative vorticity</span>
  378. <a name="l00289"></a>00289
  379. <a name="l00290"></a>00290 <span class="keyword">do</span> l = 1 , nhpp
  380. <a name="l00291"></a>00291 w = 1
  381. <a name="l00292"></a>00292 <span class="keyword">do</span> m = 1 ,ntp1
  382. <a name="l00293"></a>00293 e = w + ntp1 - m
  383. <a name="l00294"></a>00294 uds = i * dot_product(qu(w :e:2,l),pd(w: e:2))
  384. <a name="l00295"></a>00295 vds = dot_product(qv(w :e:2,l),pd(w: e:2))
  385. <a name="l00296"></a>00296 uzs = i * dot_product(qu(w :e:2,l),pz(w: e:2))
  386. <a name="l00297"></a>00297 vzs = dot_product(qv(w :e:2,l),pz(w: e:2))
  387. <a name="l00298"></a>00298 uda = i * dot_product(qu(w+1:e:2,l),pd(w+1:e:2))
  388. <a name="l00299"></a>00299 vda = dot_product(qv(w+1:e:2,l),pd(w+1:e:2))
  389. <a name="l00300"></a>00300 uza = i * dot_product(qu(w+1:e:2,l),pz(w+1:e:2))
  390. <a name="l00301"></a>00301 vza = dot_product(qv(w+1:e:2,l),pz(w+1:e:2))
  391. <a name="l00302"></a>00302 pu(m ,l) = vzs - uds + vza - uda
  392. <a name="l00303"></a>00303 pu(m+nlat,l) = -vzs - uds + vza + uda
  393. <a name="l00304"></a>00304 pv(m ,l) = -vds - uzs - vda - uza
  394. <a name="l00305"></a>00305 pv(m+nlat,l) = vds - uzs - vda + uza
  395. <a name="l00306"></a>00306 w = e + 1
  396. <a name="l00307"></a>00307 <span class="keyword">enddo</span> <span class="comment">! m</span>
  397. <a name="l00308"></a>00308 <span class="keyword">enddo</span> <span class="comment">! l</span>
  398. <a name="l00309"></a>00309 pz(2) = zsave
  399. <a name="l00310"></a>00310 return
  400. <a name="l00311"></a>00311 <span class="keyword">end</span>
  401. <a name="l00312"></a>00312
  402. <a name="l00313"></a>00313
  403. <a name="l00314"></a>00314 <span class="comment">! ================</span>
  404. <a name="l00315"></a>00315 <span class="comment">! SUBROUTINE DV2UV</span>
  405. <a name="l00316"></a>00316 <span class="comment">! ================</span>
  406. <a name="l00317"></a>00317
  407. <a name="l00318"></a><a class="code" href="legsym_8f90.html#af9cbedf7e87d9d5b2360c204237cc698">00318</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#af9cbedf7e87d9d5b2360c204237cc698">dv2uv</a>(pd,pz,pu,pv)
  408. <a name="l00319"></a>00319 use <span class="keywordflow">legsym</span>
  409. <a name="l00320"></a>00320 <span class="keyword">implicit none</span>
  410. <a name="l00321"></a>00321
  411. <a name="l00322"></a>00322 <span class="keywordtype">real</span> :: pd(2,nesp/2) <span class="comment">! Spherical harmonics of divergence</span>
  412. <a name="l00323"></a>00323 <span class="keywordtype">real</span> :: pz(2,nesp/2) <span class="comment">! Spherical harmonics of vorticity</span>
  413. <a name="l00324"></a>00324 <span class="keywordtype">real</span> :: pu(2,nlon,nhpp) <span class="comment">! Fourier coefficients of u</span>
  414. <a name="l00325"></a>00325 <span class="keywordtype">real</span> :: pv(2,nlon,nhpp) <span class="comment">! Fourier coefficients of v</span>
  415. <a name="l00326"></a>00326 <span class="keywordtype">real</span> :: zsave
  416. <a name="l00327"></a>00327 <span class="keywordtype">real</span> :: unr,uni,usr,usi,vnr,vni,vsr,vsi
  417. <a name="l00328"></a>00328 <span class="keywordtype">real</span> :: zdr,zdi,zzr,zzi
  418. <a name="l00329"></a>00329
  419. <a name="l00330"></a>00330 <span class="keywordtype">integer</span> :: l <span class="comment">! Loop index for latitude</span>
  420. <a name="l00331"></a>00331 <span class="keywordtype">integer</span> :: m <span class="comment">! Loop index for zonal wavenumber m</span>
  421. <a name="l00332"></a>00332 <span class="keywordtype">integer</span> :: n <span class="comment">! Loop index for total wavenumber n</span>
  422. <a name="l00333"></a>00333 <span class="keywordtype">integer</span> :: w <span class="comment">! Loop index for spectral mode</span>
  423. <a name="l00334"></a>00334
  424. <a name="l00335"></a>00335 pu(:,:,:) = 0.0
  425. <a name="l00336"></a>00336 pv(:,:,:) = 0.0
  426. <a name="l00337"></a>00337
  427. <a name="l00338"></a>00338 zsave = pz(1,2) <span class="comment">! Save mode(0,1) of vorticity</span>
  428. <a name="l00339"></a>00339 pz(1,2) = zsave - plavor <span class="comment">! Convert pz from absolute to relative vorticity</span>
  429. <a name="l00340"></a>00340
  430. <a name="l00341"></a>00341 <span class="keyword">do</span> l = 1 , nhpp
  431. <a name="l00342"></a>00342 w = 1
  432. <a name="l00343"></a>00343 <span class="keyword">do</span> m = 1 , ntp1
  433. <a name="l00344"></a>00344 unr = 0.0 <span class="comment">! u - north - real</span>
  434. <a name="l00345"></a>00345 uni = 0.0 <span class="comment">! u - north - imag</span>
  435. <a name="l00346"></a>00346 usr = 0.0 <span class="comment">! u - south - real</span>
  436. <a name="l00347"></a>00347 usi = 0.0 <span class="comment">! u - south - imag</span>
  437. <a name="l00348"></a>00348 vnr = 0.0 <span class="comment">! v - north - real</span>
  438. <a name="l00349"></a>00349 vni = 0.0 <span class="comment">! v - north - imag</span>
  439. <a name="l00350"></a>00350 vsr = 0.0 <span class="comment">! v - south - real</span>
  440. <a name="l00351"></a>00351 vsi = 0.0 <span class="comment">! v - south - imag</span>
  441. <a name="l00352"></a>00352
  442. <a name="l00353"></a>00353 <span class="comment">! process two modes per iteration, one symmetric (m+n = even) and one anti</span>
  443. <a name="l00354"></a>00354 <span class="comment">! we start the loop with (n=m), so the starting mode is always symmetric</span>
  444. <a name="l00355"></a>00355
  445. <a name="l00356"></a>00356 <span class="keyword">do</span> n = m , ntru , 2
  446. <a name="l00357"></a>00357 zdr = qu(w,l) * pd(1,w) <span class="comment">! symmetric mode</span>
  447. <a name="l00358"></a>00358 zdi = qu(w,l) * pd(2,w)
  448. <a name="l00359"></a>00359 zzr = qv(w,l) * pz(1,w)
  449. <a name="l00360"></a>00360 zzi = qv(w,l) * pz(2,w)
  450. <a name="l00361"></a>00361 unr = unr + zzr + zdi
  451. <a name="l00362"></a>00362 uni = uni + zzi - zdr
  452. <a name="l00363"></a>00363 usr = usr - zzr + zdi
  453. <a name="l00364"></a>00364 usi = usi - zzi - zdr
  454. <a name="l00365"></a>00365 zzr = qu(w,l) * pz(1,w)
  455. <a name="l00366"></a>00366 zzi = qu(w,l) * pz(2,w)
  456. <a name="l00367"></a>00367 zdr = qv(w,l) * pd(1,w)
  457. <a name="l00368"></a>00368 zdi = qv(w,l) * pd(2,w)
  458. <a name="l00369"></a>00369 vnr = vnr + zzi - zdr
  459. <a name="l00370"></a>00370 vni = vni - zzr - zdi
  460. <a name="l00371"></a>00371 vsr = vsr + zzi + zdr
  461. <a name="l00372"></a>00372 vsi = vsi - zzr + zdi
  462. <a name="l00373"></a>00373 w = w + 1
  463. <a name="l00374"></a>00374 zdr = qu(w,l) * pd(1,w) <span class="comment">! antisymmetric mode</span>
  464. <a name="l00375"></a>00375 zdi = qu(w,l) * pd(2,w)
  465. <a name="l00376"></a>00376 zzr = qv(w,l) * pz(1,w)
  466. <a name="l00377"></a>00377 zzi = qv(w,l) * pz(2,w)
  467. <a name="l00378"></a>00378 unr = unr + zzr + zdi
  468. <a name="l00379"></a>00379 uni = uni + zzi - zdr
  469. <a name="l00380"></a>00380 usr = usr + zzr - zdi
  470. <a name="l00381"></a>00381 usi = usi + zzi + zdr
  471. <a name="l00382"></a>00382 zzr = qu(w,l) * pz(1,w)
  472. <a name="l00383"></a>00383 zzi = qu(w,l) * pz(2,w)
  473. <a name="l00384"></a>00384 zdr = qv(w,l) * pd(1,w)
  474. <a name="l00385"></a>00385 zdi = qv(w,l) * pd(2,w)
  475. <a name="l00386"></a>00386 vnr = vnr + zzi - zdr
  476. <a name="l00387"></a>00387 vni = vni - zzr - zdi
  477. <a name="l00388"></a>00388 vsr = vsr - zzi - zdr
  478. <a name="l00389"></a>00389 vsi = vsi + zzr - zdi
  479. <a name="l00390"></a>00390 w = w + 1
  480. <a name="l00391"></a>00391 <span class="keyword">enddo</span>
  481. <a name="l00392"></a>00392 <span class="keyword">if</span> (n == ntp1) <span class="keyword">then</span> <span class="comment">! additional symmetric mode</span>
  482. <a name="l00393"></a>00393 zdr = qu(w,l) * pd(1,w) <span class="comment">! if (ntp1-m) is even</span>
  483. <a name="l00394"></a>00394 zdi = qu(w,l) * pd(2,w)
  484. <a name="l00395"></a>00395 zzr = qv(w,l) * pz(1,w)
  485. <a name="l00396"></a>00396 zzi = qv(w,l) * pz(2,w)
  486. <a name="l00397"></a>00397 unr = unr + zzr + zdi
  487. <a name="l00398"></a>00398 uni = uni + zzi - zdr
  488. <a name="l00399"></a>00399 usr = usr - zzr + zdi
  489. <a name="l00400"></a>00400 usi = usi - zzi - zdr
  490. <a name="l00401"></a>00401 zzr = qu(w,l) * pz(1,w)
  491. <a name="l00402"></a>00402 zzi = qu(w,l) * pz(2,w)
  492. <a name="l00403"></a>00403 zdr = qv(w,l) * pd(1,w)
  493. <a name="l00404"></a>00404 zdi = qv(w,l) * pd(2,w)
  494. <a name="l00405"></a>00405 vnr = vnr + zzi - zdr
  495. <a name="l00406"></a>00406 vni = vni - zzr - zdi
  496. <a name="l00407"></a>00407 vsr = vsr + zzi + zdr
  497. <a name="l00408"></a>00408 vsi = vsi - zzr + zdi
  498. <a name="l00409"></a>00409 w = w + 1
  499. <a name="l00410"></a>00410 <span class="keyword">endif</span>
  500. <a name="l00411"></a>00411 pu(1,m ,l) = unr
  501. <a name="l00412"></a>00412 pu(2,m ,l) = uni
  502. <a name="l00413"></a>00413 pu(1,m+nlat,l) = usr
  503. <a name="l00414"></a>00414 pu(2,m+nlat,l) = usi
  504. <a name="l00415"></a>00415 pv(1,m ,l) = vnr
  505. <a name="l00416"></a>00416 pv(2,m ,l) = vni
  506. <a name="l00417"></a>00417 pv(1,m+nlat,l) = vsr
  507. <a name="l00418"></a>00418 pv(2,m+nlat,l) = vsi
  508. <a name="l00419"></a>00419 <span class="keyword">enddo</span> <span class="comment">! m</span>
  509. <a name="l00420"></a>00420 <span class="keyword">enddo</span> <span class="comment">! l</span>
  510. <a name="l00421"></a>00421 pz(1,2) = zsave <span class="comment">! Restore pz to absolute vorticity</span>
  511. <a name="l00422"></a>00422 return
  512. <a name="l00423"></a>00423 <span class="keyword">end</span>
  513. <a name="l00424"></a>00424
  514. <a name="l00425"></a>00425
  515. <a name="l00426"></a>00426 <span class="comment">! =================</span>
  516. <a name="l00427"></a>00427 <span class="comment">! SUBROUTINE MKTEND</span>
  517. <a name="l00428"></a>00428 <span class="comment">! =================</span>
  518. <a name="l00429"></a>00429
  519. <a name="l00430"></a><a class="code" href="legsym_8f90.html#ab97cf272bad63e9bdd87a01317bb71c9">00430</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#ab97cf272bad63e9bdd87a01317bb71c9">mktend</a>(d,t,z,tn,fu,fv,ke,ut,vt)
  520. <a name="l00431"></a>00431 use <span class="keywordflow">legsym</span>
  521. <a name="l00432"></a>00432 <span class="keyword">implicit none</span>
  522. <a name="l00433"></a>00433
  523. <a name="l00434"></a>00434 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: tn(nlon,nhpp)
  524. <a name="l00435"></a>00435 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: fu(nlon,nhpp)
  525. <a name="l00436"></a>00436 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: fv(nlon,nhpp)
  526. <a name="l00437"></a>00437 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: ke(nlon,nhpp)
  527. <a name="l00438"></a>00438 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: ut(nlon,nhpp)
  528. <a name="l00439"></a>00439 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(in)</span> :: vt(nlon,nhpp)
  529. <a name="l00440"></a>00440
  530. <a name="l00441"></a>00441 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(out)</span> :: d(nesp/2)
  531. <a name="l00442"></a>00442 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(out)</span> :: t(nesp/2)
  532. <a name="l00443"></a>00443 <span class="keywordtype">complex</span>, <span class="keywordtype">intent(out)</span> :: z(nesp/2)
  533. <a name="l00444"></a>00444
  534. <a name="l00445"></a>00445 <span class="keywordtype">integer</span> :: l <span class="comment">! Loop index for latitude</span>
  535. <a name="l00446"></a>00446 <span class="keywordtype">integer</span> :: m <span class="comment">! Loop index for zonal wavenumber m</span>
  536. <a name="l00447"></a>00447 <span class="keywordtype">integer</span> :: w <span class="comment">! Loop index for spectral mode</span>
  537. <a name="l00448"></a>00448 <span class="keywordtype">integer</span> :: e <span class="comment">! End index for w</span>
  538. <a name="l00449"></a>00449
  539. <a name="l00450"></a>00450 <span class="keywordtype">complex</span> :: fus,fua,fvs,fva,kes,kea,tns,tna,uts,uta,vts,vta
  540. <a name="l00451"></a>00451
  541. <a name="l00452"></a>00452 d(:) = (0.0,0.0) <span class="comment">! divergence</span>
  542. <a name="l00453"></a>00453 t(:) = (0.0,0.0) <span class="comment">! temperature</span>
  543. <a name="l00454"></a>00454 z(:) = (0.0,0.0) <span class="comment">! vorticity</span>
  544. <a name="l00455"></a>00455
  545. <a name="l00456"></a>00456 <span class="keyword">do</span> l = 1 , nhpp <span class="comment">! process pairs of Nort-South latitudes</span>
  546. <a name="l00457"></a>00457 w = 1
  547. <a name="l00458"></a>00458 <span class="keyword">do</span> m = 1 , ntp1
  548. <a name="l00459"></a>00459 kes = ke(m,l) + ke(m+nlat,l) ; kea = ke(m,l) - ke(m+nlat,l)
  549. <a name="l00460"></a>00460 fvs = fv(m,l) + fv(m+nlat,l) ; fva = fv(m,l) - fv(m+nlat,l)
  550. <a name="l00461"></a>00461 fus = fu(m,l) + fu(m+nlat,l) ; fua = fu(m,l) - fu(m+nlat,l)
  551. <a name="l00462"></a>00462 uts = ut(m,l) + ut(m+nlat,l) ; uta = ut(m,l) - ut(m+nlat,l)
  552. <a name="l00463"></a>00463 vts = vt(m,l) + vt(m+nlat,l) ; vta = vt(m,l) - vt(m+nlat,l)
  553. <a name="l00464"></a>00464 tns = tn(m,l) + tn(m+nlat,l) ; tna = tn(m,l) - tn(m+nlat,l)
  554. <a name="l00465"></a>00465 e = w + ntp1 - m <span class="comment">! vector of symmetric modes</span>
  555. <a name="l00466"></a>00466 d(w:e:2) = d(w:e:2) + qq(w:e:2,l) * kes - qe(w:e:2,l) * fva + qx(w:e:2,l) * fus
  556. <a name="l00467"></a>00467 t(w:e:2) = t(w:e:2) + qe(w:e:2,l) * vta + qc(w:e:2,l) * tns - qx(w:e:2,l) * uts
  557. <a name="l00468"></a>00468 z(w:e:2) = z(w:e:2) + qe(w:e:2,l) * fua + qx(w:e:2,l) * fvs
  558. <a name="l00469"></a>00469 w = w + 1 <span class="comment">! vector of antisymmetric modes</span>
  559. <a name="l00470"></a>00470 d(w:e:2) = d(w:e:2) + qq(w:e:2,l) * kea - qe(w:e:2,l) * fvs + qx(w:e:2,l) * fua
  560. <a name="l00471"></a>00471 t(w:e:2) = t(w:e:2) + qe(w:e:2,l) * vts + qc(w:e:2,l) * tna - qx(w:e:2,l) * uta
  561. <a name="l00472"></a>00472 z(w:e:2) = z(w:e:2) + qe(w:e:2,l) * fus + qx(w:e:2,l) * fva
  562. <a name="l00473"></a>00473 w = e + 1
  563. <a name="l00474"></a>00474 <span class="keyword">enddo</span> <span class="comment">! m</span>
  564. <a name="l00475"></a>00475 <span class="keyword">enddo</span> <span class="comment">! l</span>
  565. <a name="l00476"></a>00476 return
  566. <a name="l00477"></a>00477 <span class="keyword">end</span>
  567. <a name="l00478"></a>00478
  568. <a name="l00479"></a>00479
  569. <a name="l00480"></a>00480 <span class="comment">! ==================</span>
  570. <a name="l00481"></a>00481 <span class="comment">! SUBROUTINE REG2ALT</span>
  571. <a name="l00482"></a>00482 <span class="comment">! ==================</span>
  572. <a name="l00483"></a>00483
  573. <a name="l00484"></a><a class="code" href="legsym_8f90.html#a4a468562c0549b4ca3ec6ea34f87545a">00484</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a4a468562c0549b4ca3ec6ea34f87545a">reg2alt</a>(pr,klev)
  574. <a name="l00485"></a>00485 use <span class="keywordflow">legsym</span>
  575. <a name="l00486"></a>00486 <span class="keyword">implicit none</span>
  576. <a name="l00487"></a>00487
  577. <a name="l00488"></a>00488 <span class="keywordtype">real</span> :: pr(nlon,nlat,klev)
  578. <a name="l00489"></a>00489 <span class="keywordtype">real</span> :: pa(nlon,nlat,klev)
  579. <a name="l00490"></a>00490
  580. <a name="l00491"></a>00491 <span class="keywordtype">integer</span> :: jlat
  581. <a name="l00492"></a>00492 <span class="keywordtype">integer</span> :: klev
  582. <a name="l00493"></a>00493
  583. <a name="l00494"></a>00494 <span class="keyword">do</span> jlat = 1 , nlat / 2
  584. <a name="l00495"></a>00495 pa(:,2*jlat-1,:) = pr(:,jlat ,:)
  585. <a name="l00496"></a>00496 pa(:,2*jlat ,:) = pr(:,nlat-jlat+1,:)
  586. <a name="l00497"></a>00497 <span class="keyword">enddo</span>
  587. <a name="l00498"></a>00498
  588. <a name="l00499"></a>00499 pr = pa
  589. <a name="l00500"></a>00500
  590. <a name="l00501"></a>00501 return
  591. <a name="l00502"></a>00502 <span class="keyword">end</span>
  592. <a name="l00503"></a>00503
  593. <a name="l00504"></a>00504
  594. <a name="l00505"></a>00505 <span class="comment">! ==================</span>
  595. <a name="l00506"></a>00506 <span class="comment">! SUBROUTINE ALT2REG</span>
  596. <a name="l00507"></a>00507 <span class="comment">! ==================</span>
  597. <a name="l00508"></a>00508
  598. <a name="l00509"></a><a class="code" href="legsym_8f90.html#a308819246e409c8dbe1e778d304ef415">00509</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a308819246e409c8dbe1e778d304ef415">alt2reg</a>(pa,klev)
  599. <a name="l00510"></a>00510 use <span class="keywordflow">legsym</span>
  600. <a name="l00511"></a>00511 <span class="keyword">implicit none</span>
  601. <a name="l00512"></a>00512
  602. <a name="l00513"></a>00513 <span class="keywordtype">real</span> :: pa(nlon,nlat,klev)
  603. <a name="l00514"></a>00514 <span class="keywordtype">real</span> :: pr(nlon,nlat,klev)
  604. <a name="l00515"></a>00515
  605. <a name="l00516"></a>00516 <span class="keywordtype">integer</span> :: jlat
  606. <a name="l00517"></a>00517 <span class="keywordtype">integer</span> :: klev
  607. <a name="l00518"></a>00518
  608. <a name="l00519"></a>00519 <span class="keyword">do</span> jlat = 1 , nlat / 2
  609. <a name="l00520"></a>00520 pr(:,jlat ,:) = pa(:,2*jlat-1,:)
  610. <a name="l00521"></a>00521 pr(:,nlat-jlat+1,:) = pa(:,2*jlat ,:)
  611. <a name="l00522"></a>00522 <span class="keyword">enddo</span>
  612. <a name="l00523"></a>00523
  613. <a name="l00524"></a>00524 pa = pr
  614. <a name="l00525"></a>00525
  615. <a name="l00526"></a>00526 return
  616. <a name="l00527"></a>00527 <span class="keyword">end</span>
  617. <a name="l00528"></a>00528
  618. <a name="l00529"></a>00529
  619. <a name="l00530"></a>00530 <span class="comment">! ================</span>
  620. <a name="l00531"></a>00531 <span class="comment">! SUBROUTINE ALTCS</span>
  621. <a name="l00532"></a>00532 <span class="comment">! ================</span>
  622. <a name="l00533"></a>00533
  623. <a name="l00534"></a><a class="code" href="legsym_8f90.html#a6ba5b0b99819bcbad73f2e2eb49c62bb">00534</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#a6ba5b0b99819bcbad73f2e2eb49c62bb">altcs</a>(pcs)
  624. <a name="l00535"></a>00535 use <span class="keywordflow">legsym</span>
  625. <a name="l00536"></a>00536 <span class="keyword">implicit none</span>
  626. <a name="l00537"></a>00537 <span class="keywordtype">real</span> :: pcs(nlat,nlev)
  627. <a name="l00538"></a>00538 <span class="keywordtype">real</span> :: pal(nlat,nlev)
  628. <a name="l00539"></a>00539 <span class="keywordtype">integer</span> :: jlat
  629. <a name="l00540"></a>00540
  630. <a name="l00541"></a>00541 <span class="keyword">do</span> jlat = 1 , nlat / 2
  631. <a name="l00542"></a>00542 pal(jlat ,:) = pcs(2*jlat-1,:)
  632. <a name="l00543"></a>00543 pal(nlat-jlat+1,:) = pcs(2*jlat ,:)
  633. <a name="l00544"></a>00544 <span class="keyword">enddo</span>
  634. <a name="l00545"></a>00545
  635. <a name="l00546"></a>00546 pcs = pal
  636. <a name="l00547"></a>00547
  637. <a name="l00548"></a>00548 return
  638. <a name="l00549"></a>00549 <span class="keyword">end</span>
  639. <a name="l00550"></a>00550
  640. <a name="l00551"></a>00551
  641. <a name="l00552"></a>00552 <span class="comment">! =================</span>
  642. <a name="l00553"></a>00553 <span class="comment">! SUBROUTINE ALTLAT </span>
  643. <a name="l00554"></a>00554 <span class="comment">! =================</span>
  644. <a name="l00555"></a>00555
  645. <a name="l00556"></a><a class="code" href="legsym_8f90.html#ae810767bcafdac840ab48c420efcb49a">00556</a> <span class="keyword">subroutine </span><a class="code" href="legsym_8f90.html#ae810767bcafdac840ab48c420efcb49a">altlat</a>(pr,klat)
  646. <a name="l00557"></a>00557 <span class="keyword">implicit none</span>
  647. <a name="l00558"></a>00558 <span class="keywordtype">integer</span> :: jlat
  648. <a name="l00559"></a>00559 <span class="keywordtype">integer</span> :: klat
  649. <a name="l00560"></a>00560 <span class="keywordtype">real</span> :: pr(klat) <span class="comment">! regular grid</span>
  650. <a name="l00561"></a>00561 <span class="keywordtype">real</span> :: pa(klat) <span class="comment">! alternating grid</span>
  651. <a name="l00562"></a>00562
  652. <a name="l00563"></a>00563 <span class="keyword">do</span> jlat = 1 , klat / 2
  653. <a name="l00564"></a>00564 pa(2*jlat-1) = pr(jlat )
  654. <a name="l00565"></a>00565 pa(2*jlat ) = pr(klat-jlat+1)
  655. <a name="l00566"></a>00566 <span class="keyword">enddo</span>
  656. <a name="l00567"></a>00567
  657. <a name="l00568"></a>00568 pr(:) = pa(:)
  658. <a name="l00569"></a>00569
  659. <a name="l00570"></a>00570 return
  660. <a name="l00571"></a>00571 <span class="keyword">end</span>
  661. </pre></div></div>
  662. </div>
  663. <div id="nav-path" class="navpath">
  664. <ul>
  665. <li class="navelem"><a class="el" href="legsym_8f90.html">legsym.f90</a> </li>
  666. <!-- window showing the filter options -->
  667. <div id="MSearchSelectWindow"
  668. onmouseover="return searchBox.OnSearchSelectShow()"
  669. onmouseout="return searchBox.OnSearchSelectHide()"
  670. onkeydown="return searchBox.OnSearchSelectKey(event)">
  671. <a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Defines</a></div>
  672. <!-- iframe showing the search results (closed by default) -->
  673. <div id="MSearchResultsWindow">
  674. <iframe src="javascript:void(0)" frameborder="0"
  675. name="MSearchResults" id="MSearchResults">
  676. </iframe>
  677. </div>
  678. <li class="footer">Generated on Wed Sep 21 2011 12:35:46 for PUMA by
  679. <a href="http://www.doxygen.org/index.html">
  680. <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.5.1 </li>
  681. </ul>
  682. </div>
  683. </body>
  684. </html>