123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517 |
- <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
- <html xmlns="http://www.w3.org/1999/xhtml">
- <head>
- <meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
- <title>PUMA: /Users/home/WC/puma/src/mpimod_multi.f90 Source File</title>
- <link href="tabs.css" rel="stylesheet" type="text/css"/>
- <link href="doxygen.css" rel="stylesheet" type="text/css" />
- <link href="navtree.css" rel="stylesheet" type="text/css"/>
- <script type="text/javascript" src="jquery.js"></script>
- <script type="text/javascript" src="resize.js"></script>
- <script type="text/javascript" src="navtree.js"></script>
- <script type="text/javascript">
- $(document).ready(initResizable);
- </script>
- <link href="search/search.css" rel="stylesheet" type="text/css"/>
- <script type="text/javascript" src="search/search.js"></script>
- <script type="text/javascript">
- $(document).ready(function() { searchBox.OnSelectItem(0); });
- </script>
- </head>
- <body>
- <div id="top"><!-- do not remove this div! -->
- <div id="titlearea">
- <table cellspacing="0" cellpadding="0">
- <tbody>
- <tr style="height: 56px;">
-
- <td id="projectlogo"><img alt="Logo" src="puma103.jpg"/></td>
-
-
- <td style="padding-left: 0.5em;">
- <div id="projectname">PUMA
-  <span id="projectnumber">219</span>
- </div>
- <div id="projectbrief">Portable University Model of the Atmosphere</div>
- </td>
-
-
-
- </tr>
- </tbody>
- </table>
- </div>
- <!-- Generated by Doxygen 1.7.5.1 -->
- <script type="text/javascript">
- var searchBox = new SearchBox("searchBox", "search",false,'Search');
- </script>
- <div id="navrow1" class="tabs">
- <ul class="tablist">
- <li><a href="index.html"><span>Main Page</span></a></li>
- <li><a href="annotated.html"><span>Data Types List</span></a></li>
- <li class="current"><a href="files.html"><span>Files</span></a></li>
- <li>
- <div id="MSearchBox" class="MSearchBoxInactive">
- <span class="left">
- <img id="MSearchSelect" src="search/mag_sel.png"
- onmouseover="return searchBox.OnSearchSelectShow()"
- onmouseout="return searchBox.OnSearchSelectHide()"
- alt=""/>
- <input type="text" id="MSearchField" value="Search" accesskey="S"
- onfocus="searchBox.OnSearchFieldFocus(true)"
- onblur="searchBox.OnSearchFieldFocus(false)"
- onkeyup="searchBox.OnSearchFieldChange(event)"/>
- </span><span class="right">
- <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
- </span>
- </div>
- </li>
- </ul>
- </div>
- <div id="navrow2" class="tabs2">
- <ul class="tablist">
- <li><a href="files.html"><span>File List</span></a></li>
- <li><a href="globals.html"><span>File Members</span></a></li>
- </ul>
- </div>
- </div>
- <div id="side-nav" class="ui-resizable side-nav-resizable">
- <div id="nav-tree">
- <div id="nav-tree-contents">
- </div>
- </div>
- <div id="splitbar" style="-moz-user-select:none;"
- class="ui-resizable-handle">
- </div>
- </div>
- <script type="text/javascript">
- initNavTree('mpimod__multi_8f90.html','');
- </script>
- <div id="doc-content">
- <div class="header">
- <div class="headertitle">
- <div class="title">/Users/home/WC/puma/src/mpimod_multi.f90</div> </div>
- </div>
- <div class="contents">
- <a href="mpimod__multi_8f90.html">Go to the documentation of this file.</a><div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 <span class="keyword">module</span> <a class="code" href="classmpimod.html">mpimod</a>
- <a name="l00002"></a>00002 use <span class="keywordflow">pumamod</span>
- <a name="l00003"></a>00003 use <span class="keywordflow">mpi</span>
- <a name="l00004"></a>00004
- <a name="l00005"></a>00005 <span class="keywordtype">integer</span> :: mpi_itype = MPI_INTEGER4
- <a name="l00006"></a>00006 <span class="keywordtype">integer</span> :: mpi_rtype = MPI_REAL4
- <a name="l00007"></a>00007 <span class="keywordtype">integer</span> :: mpi_ltype = MPI_LOGICAL
- <a name="l00008"></a>00008
- <a name="l00009"></a>00009 <span class="keyword"> end module mpimod</span>
- <a name="l00010"></a>00010
- <a name="l00011"></a>00011
- <a name="l00012"></a>00012 <span class="comment">! ==================</span>
- <a name="l00013"></a>00013 <span class="comment">! SUBROUTINE MPSTART</span>
- <a name="l00014"></a>00014 <span class="comment">! ==================</span>
- <a name="l00015"></a>00015
- <a name="l00016"></a><a class="code" href="mpimod__multi_8f90.html#a41bbd9334a3d0412c73399d699bbb237">00016</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a41bbd9334a3d0412c73399d699bbb237">mpstart</a> <span class="comment">! initialization</span>
- <a name="l00017"></a>00017 use <span class="keywordflow">mpimod</span>
- <a name="l00018"></a>00018 <span class="keywordtype">integer</span> :: itest = 0
- <a name="l00019"></a>00019 <span class="keywordtype">real</span> :: rtest = 0.0
- <a name="l00020"></a>00020
- <a name="l00021"></a>00021 <span class="keyword">if</span> (kind(itest) == 8) mpi_itype = MPI_INTEGER8
- <a name="l00022"></a>00022 <span class="keyword">if</span> (kind(rtest) == 8) mpi_rtype = MPI_REAL8
- <a name="l00023"></a>00023
- <a name="l00024"></a>00024 call mpi_init(mrinfo)
- <a name="l00025"></a>00025 mrworld=MPI_COMM_WORLD
- <a name="l00026"></a>00026
- <a name="l00027"></a>00027 call mpi_comm_size(mrworld,mrnum,mrinfo)
- <a name="l00028"></a>00028 call mpi_comm_rank(mrworld,mrpid,mrinfo)
- <a name="l00029"></a>00029 <span class="keyword">allocate</span>(ympname(mrnum)) ; ympname(:) = <span class="stringliteral">' '</span> <span class="comment">! process names</span>
- <a name="l00030"></a>00030 call mpi_get_processor_name(ympname(1),ilen,mrinfo)
- <a name="l00031"></a>00031
- <a name="l00032"></a>00032 call mpi_gather(ympname,80,mpi_character, &
- <a name="l00033"></a>00033 ympname,80,mpi_character, &
- <a name="l00034"></a>00034 nroot,mrworld,mrinfo)
- <a name="l00035"></a>00035 return
- <a name="l00036"></a>00036 <span class="keyword"> end subroutine mpstart</span>
- <a name="l00037"></a>00037
- <a name="l00038"></a>00038
- <a name="l00039"></a>00039 <span class="comment">! =======================</span>
- <a name="l00040"></a>00040 <span class="comment">! subroutine mrdimensions</span>
- <a name="l00041"></a>00041 <span class="comment">! =======================</span>
- <a name="l00042"></a>00042
- <a name="l00043"></a><a class="code" href="mpimod__multi_8f90.html#acb4a2403b5f65a70e7e5ff01ea4577f7">00043</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#acb4a2403b5f65a70e7e5ff01ea4577f7">mrdimensions</a>
- <a name="l00044"></a>00044 use <span class="keywordflow">mpimod</span>
- <a name="l00045"></a>00045
- <a name="l00046"></a>00046 <span class="keyword">allocate</span>(mrtru(mrnum)) <span class="comment">! all truncations</span>
- <a name="l00047"></a>00047 mrtru(1) = ntru
- <a name="l00048"></a>00048 call mpi_gather(mrtru,1,mpi_itype, &
- <a name="l00049"></a>00049 mrtru,1,mpi_itype, &
- <a name="l00050"></a>00050 nroot,mrworld,mrinfo)
- <a name="l00051"></a>00051
- <a name="l00052"></a>00052 call <a class="code" href="mpimod__multi_8f90.html#a7fd328bc0c83a30f521615978d5d44b2">mrbcin</a>(mrtru,mrnum) <span class="comment">! broadcast truncations</span>
- <a name="l00053"></a>00053 mintru = minval(mrtru(1:mrnum))
- <a name="l00054"></a>00054 mrdim = (mintru+1) * (mintru+2)
- <a name="l00055"></a>00055 return
- <a name="l00056"></a>00056 <span class="keyword"> end</span>
- <a name="l00057"></a>00057
- <a name="l00058"></a>00058 <span class="comment">! =================</span>
- <a name="l00059"></a>00059 <span class="comment">! subroutine mpstop</span>
- <a name="l00060"></a>00060 <span class="comment">! =================</span>
- <a name="l00061"></a>00061
- <a name="l00062"></a><a class="code" href="mpimod__multi_8f90.html#ac80e83b9bc0a4b459fed5f3b79cfafa0">00062</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#ac80e83b9bc0a4b459fed5f3b79cfafa0">mpstop</a>
- <a name="l00063"></a>00063 use <span class="keywordflow">mpimod</span>
- <a name="l00064"></a>00064
- <a name="l00065"></a>00065 call mpi_barrier(mrworld,mrinfo)
- <a name="l00066"></a>00066 call mpi_finalize(mrinfo)
- <a name="l00067"></a>00067
- <a name="l00068"></a>00068 return
- <a name="l00069"></a>00069 <span class="keyword"> end subroutine mpstop</span>
- <a name="l00070"></a>00070
- <a name="l00071"></a>00071
- <a name="l00072"></a>00072 <span class="comment">! ================</span>
- <a name="l00073"></a>00073 <span class="comment">! subroutine mrsum</span>
- <a name="l00074"></a>00074 <span class="comment">! ================</span>
- <a name="l00075"></a>00075
- <a name="l00076"></a><a class="code" href="mpimod__multi_8f90.html#a5d2bb9cfe68e5feb6de6b359f04398e3">00076</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a5d2bb9cfe68e5feb6de6b359f04398e3">mrsum</a>(k)
- <a name="l00077"></a>00077 use <span class="keywordflow">mpimod</span>
- <a name="l00078"></a>00078 <span class="keywordtype">integer</span> :: k,n
- <a name="l00079"></a>00079 <span class="comment">! print *,'mrsum b[',mrpid,']',k</span>
- <a name="l00080"></a>00080 call mpi_allreduce(k,n,1,mpi_itype,MPI_SUM,mrworld,mrinfo)
- <a name="l00081"></a>00081 <span class="comment">! print *,'mrsum a[',mrpid,']',n</span>
- <a name="l00082"></a>00082 k = n
- <a name="l00083"></a>00083 return
- <a name="l00084"></a>00084 <span class="keyword"> end</span>
- <a name="l00085"></a>00085
- <a name="l00086"></a>00086
- <a name="l00087"></a>00087 <span class="comment">! =================</span>
- <a name="l00088"></a>00088 <span class="comment">! subroutine mrdiff</span>
- <a name="l00089"></a>00089 <span class="comment">! =================</span>
- <a name="l00090"></a>00090
- <a name="l00091"></a><a class="code" href="mpimod__multi_8f90.html#a81baf40d09fd23640c01877173e0813a">00091</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#af3212261e3ce775f26d09859c337b760">mrdiff</a>(pa,pd,kesp,klev)
- <a name="l00092"></a>00092 use <span class="keywordflow">mpimod</span>
- <a name="l00093"></a>00093
- <a name="l00094"></a>00094 <span class="keywordtype">real</span> :: pa(kesp,klev)
- <a name="l00095"></a>00095 <span class="keywordtype">real</span> :: pd(kesp,klev)
- <a name="l00096"></a>00096 <span class="keywordtype">real</span> :: px(mrdim)
- <a name="l00097"></a>00097 <span class="keywordtype">real</span> :: py(mrdim)
- <a name="l00098"></a>00098
- <a name="l00099"></a>00099 <span class="keyword">do</span> jlev = 1 , klev
- <a name="l00100"></a>00100 <span class="keyword">if</span> (kesp == mrdim) <span class="keyword">then</span>
- <a name="l00101"></a>00101 px(:) = pa(:,jlev)
- <a name="l00102"></a>00102 <span class="keyword">else</span>
- <a name="l00103"></a>00103 call <a class="code" href="mpimod__multi_8f90.html#a63020ca5fd4738534a8ad6f1ce97133d">mrtrunc</a>(pa(1,jlev),ntru,px,mintru)
- <a name="l00104"></a>00104 <span class="keyword">endif</span>
- <a name="l00105"></a>00105 call mpi_allreduce(px,py,mrdim,mpi_rtype,MPI_SUM,mrworld,mrinfo)
- <a name="l00106"></a>00106 py(:) = py(:) - 2.0 * px(:)
- <a name="l00107"></a>00107 <span class="keyword">if</span> (kesp == mrdim) <span class="keyword">then</span>
- <a name="l00108"></a>00108 pd(:,jlev) = py(:)
- <a name="l00109"></a>00109 <span class="keyword">else</span>
- <a name="l00110"></a>00110 call <a class="code" href="mpimod__multi_8f90.html#a7878a87944a3261177c0321cb0d01174">mrexpand</a>(py,mintru,pd(1,jlev),ntru)
- <a name="l00111"></a>00111 <span class="keyword">endif</span>
- <a name="l00112"></a>00112 <span class="keyword">enddo</span>
- <a name="l00113"></a>00113 return
- <a name="l00114"></a>00114 <span class="keyword"> end</span>
- <a name="l00115"></a>00115
- <a name="l00116"></a>00116
- <a name="l00117"></a>00117 <span class="comment">! =================</span>
- <a name="l00118"></a>00118 <span class="comment">! subroutine mrbcin</span>
- <a name="l00119"></a>00119 <span class="comment">! =================</span>
- <a name="l00120"></a>00120
- <a name="l00121"></a><a class="code" href="mpimod__multi_8f90.html#a7fd328bc0c83a30f521615978d5d44b2">00121</a> <span class="keyword">subroutine </span><a class="code" href="mpimod__multi_8f90.html#a7fd328bc0c83a30f521615978d5d44b2">mrbcin</a>(k,n)
- <a name="l00122"></a>00122 use <span class="keywordflow">mpimod</span>
- <a name="l00123"></a>00123 <span class="keywordtype">integer</span> :: k(n)
- <a name="l00124"></a>00124
- <a name="l00125"></a>00125 call mpi_bcast(k,n,mpi_itype,NROOT,mrworld,mrinfo)
- <a name="l00126"></a>00126
- <a name="l00127"></a>00127 return
- <a name="l00128"></a>00128 <span class="keyword"> end</span>
- <a name="l00129"></a>00129
- <a name="l00130"></a>00130
- <a name="l00131"></a>00131 <span class="comment">! ==================</span>
- <a name="l00132"></a>00132 <span class="comment">! subroutine mrtrunc</span>
- <a name="l00133"></a>00133 <span class="comment">! ==================</span>
- <a name="l00134"></a>00134
- <a name="l00135"></a>00135 <span class="comment">! Reduce truncation from nthi to ntlo</span>
- <a name="l00136"></a>00136
- <a name="l00137"></a><a class="code" href="mpimod__multi_8f90.html#a63020ca5fd4738534a8ad6f1ce97133d">00137</a> <span class="keyword">subroutine </span><a class="code" href="mpimod__multi_8f90.html#a63020ca5fd4738534a8ad6f1ce97133d">mrtrunc</a>(sphi,nthi,splo,ntlo)
- <a name="l00138"></a>00138 <span class="keywordtype">complex</span> :: sphi(*)
- <a name="l00139"></a>00139 <span class="keywordtype">complex</span> :: splo(*)
- <a name="l00140"></a>00140
- <a name="l00141"></a>00141 jl = 1
- <a name="l00142"></a>00142 jh = 1
- <a name="l00143"></a>00143 <span class="keyword">do</span> jm = 0 , ntlo
- <a name="l00144"></a>00144 jn = ntlo - jm
- <a name="l00145"></a>00145 splo(jl:jl+jn) = sphi(jh:jh+jn)
- <a name="l00146"></a>00146 jl = jl + ntlo - jm + 1
- <a name="l00147"></a>00147 jh = jh + nthi - jm + 1
- <a name="l00148"></a>00148 <span class="keyword">enddo</span>
- <a name="l00149"></a>00149 return
- <a name="l00150"></a>00150 <span class="keyword"> end</span>
- <a name="l00151"></a>00151
- <a name="l00152"></a>00152
- <a name="l00153"></a>00153 <span class="comment">! ===================</span>
- <a name="l00154"></a>00154 <span class="comment">! subroutine mrexpand</span>
- <a name="l00155"></a>00155 <span class="comment">! ===================</span>
- <a name="l00156"></a>00156
- <a name="l00157"></a>00157 <span class="comment">! Expand truncation from ntlo to nthi</span>
- <a name="l00158"></a>00158
- <a name="l00159"></a><a class="code" href="mpimod__multi_8f90.html#a7878a87944a3261177c0321cb0d01174">00159</a> <span class="keyword">subroutine </span><a class="code" href="mpimod__multi_8f90.html#a7878a87944a3261177c0321cb0d01174">mrexpand</a>(splo,ntlo,sphi,nthi)
- <a name="l00160"></a>00160 <span class="keywordtype">complex</span> :: sphi(*)
- <a name="l00161"></a>00161 <span class="keywordtype">complex</span> :: splo(*)
- <a name="l00162"></a>00162
- <a name="l00163"></a>00163 sphi(1:nthi) = (0.0,0.0)
- <a name="l00164"></a>00164 jl = 1
- <a name="l00165"></a>00165 jh = 1
- <a name="l00166"></a>00166 <span class="keyword">do</span> jm = 0 , ntlo
- <a name="l00167"></a>00167 jn = ntlo - jm
- <a name="l00168"></a>00168 sphi(jh:jh+jn) = splo(jl:jl+jn)
- <a name="l00169"></a>00169 jl = jl + ntlo - jm + 1
- <a name="l00170"></a>00170 jh = jh + nthi - jm + 1
- <a name="l00171"></a>00171 <span class="keyword">enddo</span>
- <a name="l00172"></a>00172 return
- <a name="l00173"></a>00173 <span class="keyword"> end</span>
- <a name="l00174"></a>00174
- <a name="l00175"></a>00175
- <a name="l00176"></a><a class="code" href="mpimod__multi_8f90.html#a89982355acc98319bfc191dab28da805">00176</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a89982355acc98319bfc191dab28da805">mpbci</a>(k) ! broadcast 1 integer
- <a name="l00177"></a>00177 return
- <a name="l00178"></a>00178 <span class="keyword"> end</span>
- <a name="l00179"></a>00179
- <a name="l00180"></a><a class="code" href="mpimod__multi_8f90.html#a85cfae5acde5c37604edf690e9c2f7cf">00180</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a85cfae5acde5c37604edf690e9c2f7cf">mpbcin</a>(k,n) ! broadcast n integer
- <a name="l00181"></a>00181 <span class="keywordtype">integer</span> :: k(n)
- <a name="l00182"></a>00182 return
- <a name="l00183"></a>00183 <span class="keyword"> end</span>
- <a name="l00184"></a>00184
- <a name="l00185"></a><a class="code" href="mpimod__multi_8f90.html#aded092db7f8071a727e2e96887702ca7">00185</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#aded092db7f8071a727e2e96887702ca7">mpbcr</a>(p) ! broadcast 1 real
- <a name="l00186"></a>00186 return
- <a name="l00187"></a>00187 <span class="keyword"> end</span>
- <a name="l00188"></a>00188
- <a name="l00189"></a><a class="code" href="mpimod__multi_8f90.html#af2a0a009162180d4abb1daa1bad60cf2">00189</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#af2a0a009162180d4abb1daa1bad60cf2">mpbcrn</a>(p,n) ! broadcast n real
- <a name="l00190"></a>00190 <span class="keywordtype">real</span> :: p(n)
- <a name="l00191"></a>00191 return
- <a name="l00192"></a>00192 <span class="keyword"> end</span>
- <a name="l00193"></a>00193
- <a name="l00194"></a><a class="code" href="mpimod__multi_8f90.html#a40b910e38273e7f3c9dc4ed36d3e67a0">00194</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a40b910e38273e7f3c9dc4ed36d3e67a0">mpbcl</a>(k) ! broadcast 1 logical
- <a name="l00195"></a>00195 <span class="keywordtype">logical</span> :: k
- <a name="l00196"></a>00196 return
- <a name="l00197"></a>00197 <span class="keyword"> end</span>
- <a name="l00198"></a>00198
- <a name="l00199"></a><a class="code" href="mpimod__multi_8f90.html#a8338d8609afcefbb1faa41f353c10ef9">00199</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a8338d8609afcefbb1faa41f353c10ef9">mpscin</a>(k,n) ! scatter n integer
- <a name="l00200"></a>00200 <span class="keywordtype">integer</span> :: k(n)
- <a name="l00201"></a>00201 return
- <a name="l00202"></a>00202 <span class="keyword"> end</span>
- <a name="l00203"></a>00203
- <a name="l00204"></a><a class="code" href="mpimod__multi_8f90.html#a1504cf64a1ffc198a8a1fe54ba00d775">00204</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a1504cf64a1ffc198a8a1fe54ba00d775">mpscrn</a>(p,n) ! scatter n real
- <a name="l00205"></a>00205 <span class="keywordtype">real</span> :: p(n)
- <a name="l00206"></a>00206 return
- <a name="l00207"></a>00207 <span class="keyword"> end</span>
- <a name="l00208"></a>00208
- <a name="l00209"></a><a class="code" href="mpimod__multi_8f90.html#a3d2a5d231fd9527bcbc1fde327326922">00209</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a3d2a5d231fd9527bcbc1fde327326922">mpscdn</a>(p,n) ! scatter n double precision
- <a name="l00210"></a>00210 <span class="keywordtype">real (kind=8)</span> :: p(n)
- <a name="l00211"></a>00211 return
- <a name="l00212"></a>00212 <span class="keyword"> end</span>
- <a name="l00213"></a>00213
- <a name="l00214"></a><a class="code" href="mpimod__multi_8f90.html#a0c5adf4e8c7e39cf5a1038a1d34ebf30">00214</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a0c5adf4e8c7e39cf5a1038a1d34ebf30">mpscsp</a>(pf,pp,klev) ! scatter spectral fields
- <a name="l00215"></a>00215 use <span class="keywordflow">pumamod</span>
- <a name="l00216"></a>00216 <span class="keywordtype">real</span> pf(nesp,klev)
- <a name="l00217"></a>00217 <span class="keywordtype">real</span> pp(nspp,klev)
- <a name="l00218"></a>00218 pp(1:nspp,1:klev) = pf(1:nspp,1:klev)
- <a name="l00219"></a>00219 return
- <a name="l00220"></a>00220 <span class="keyword"> end</span>
- <a name="l00221"></a>00221
- <a name="l00222"></a><a class="code" href="mpimod__multi_8f90.html#ac66e76c6144dfeadbc03bc5817553250">00222</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#ac66e76c6144dfeadbc03bc5817553250">mpscgp</a>(pf,pp,klev) ! scatter gridpoint fields
- <a name="l00223"></a>00223 use <span class="keywordflow">pumamod</span>
- <a name="l00224"></a>00224 <span class="keywordtype">real</span> pf(nlon*nlat,klev)
- <a name="l00225"></a>00225 <span class="keywordtype">real</span> pp(nhor,klev)
- <a name="l00226"></a>00226 pp(1:nhor,1:klev) = pf(1:nhor,1:klev)
- <a name="l00227"></a>00227 return
- <a name="l00228"></a>00228 <span class="keyword"> end</span>
- <a name="l00229"></a>00229
- <a name="l00230"></a><a class="code" href="mpimod__multi_8f90.html#ac053a575b1230f8e4a296164dba5ab27">00230</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#ac053a575b1230f8e4a296164dba5ab27">mpgasp</a>(pf,pp,klev) ! gather spectral fields
- <a name="l00231"></a>00231 use <span class="keywordflow">pumamod</span>
- <a name="l00232"></a>00232 <span class="keywordtype">real</span> pf(nesp,klev)
- <a name="l00233"></a>00233 <span class="keywordtype">real</span> pp(nspp,klev)
- <a name="l00234"></a>00234 pf(1:nspp,1:klev) = pp(1:nspp,1:klev)
- <a name="l00235"></a>00235 return
- <a name="l00236"></a>00236 <span class="keyword"> end</span>
- <a name="l00237"></a>00237
- <a name="l00238"></a><a class="code" href="mpimod__multi_8f90.html#aaa1210298789f4fd7b7702c276eb80a9">00238</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#aaa1210298789f4fd7b7702c276eb80a9">mpgagp</a>(pf,pp,klev) ! gather gridpoint fields
- <a name="l00239"></a>00239 use <span class="keywordflow">pumamod</span>
- <a name="l00240"></a>00240 <span class="keywordtype">real</span> pf(nhor,klev)
- <a name="l00241"></a>00241 <span class="keywordtype">real</span> pp(nhor,klev)
- <a name="l00242"></a>00242 pf = pp
- <a name="l00243"></a>00243 return
- <a name="l00244"></a>00244 <span class="keyword"> end</span>
- <a name="l00245"></a>00245
- <a name="l00246"></a><a class="code" href="mpimod__multi_8f90.html#a5aef7e33503e0c46b1d8c0b984c398d1">00246</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a5aef7e33503e0c46b1d8c0b984c398d1">mpgacs</a>(pcs) ! gather cross sections
- <a name="l00247"></a>00247 return
- <a name="l00248"></a>00248 <span class="keyword"> end</span>
- <a name="l00249"></a>00249
- <a name="l00250"></a><a class="code" href="mpimod__multi_8f90.html#a54cf45feb57177de8eaab2e6b01a7aa2">00250</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a54cf45feb57177de8eaab2e6b01a7aa2">mpgallsp</a>(pf,pp,klev) ! gather spectral to all
- <a name="l00251"></a>00251 use <span class="keywordflow">pumamod</span>
- <a name="l00252"></a>00252 <span class="keywordtype">real</span> pf(nesp,klev)
- <a name="l00253"></a>00253 <span class="keywordtype">real</span> pp(nspp,klev)
- <a name="l00254"></a>00254 pf(1:nspp,1:klev) = pp(1:nspp,1:klev)
- <a name="l00255"></a>00255 return
- <a name="l00256"></a>00256 <span class="keyword"> end</span>
- <a name="l00257"></a>00257
- <a name="l00258"></a><a class="code" href="mpimod__multi_8f90.html#af894efd9525c935f22415e017dcbc482">00258</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#af894efd9525c935f22415e017dcbc482">mpsum</a>(psp,klev) ! sum spectral fields
- <a name="l00259"></a>00259 return
- <a name="l00260"></a>00260 <span class="keyword"> end</span>
- <a name="l00261"></a>00261
- <a name="l00262"></a><a class="code" href="mpimod__multi_8f90.html#a75a681a8d4b9ab5ba0d4fa97f909647b">00262</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a75a681a8d4b9ab5ba0d4fa97f909647b">mpsumsc</a>(psf,psp,klev) ! sum & scatter spectral
- <a name="l00263"></a>00263 use <span class="keywordflow">pumamod</span>
- <a name="l00264"></a>00264 <span class="keywordtype">real</span> psf(nesp,klev)
- <a name="l00265"></a>00265 <span class="keywordtype">real</span> psp(nspp,klev)
- <a name="l00266"></a>00266 psp(1:nspp,1:klev) = psf(1:nspp,1:klev)
- <a name="l00267"></a>00267 return
- <a name="l00268"></a>00268 <span class="keyword"> end</span>
- <a name="l00269"></a>00269
- <a name="l00270"></a><a class="code" href="mpimod__multi_8f90.html#af2111ef6d5b772479a74e94d351440f0">00270</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#af2111ef6d5b772479a74e94d351440f0">mpsumr</a>(pr,kdim) ! sum kdim reals
- <a name="l00271"></a>00271 return
- <a name="l00272"></a>00272 <span class="keyword"> end subroutine mpsumr</span>
- <a name="l00273"></a>00273
- <a name="l00274"></a><a class="code" href="mpimod__multi_8f90.html#ad703e6ecd123e9b8280322e402d57d20">00274</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#ad703e6ecd123e9b8280322e402d57d20">mpsumbcr</a>(pr,kdim) ! sum & broadcast kdim reals
- <a name="l00275"></a>00275 return
- <a name="l00276"></a>00276 <span class="keyword"> end</span>
- <a name="l00277"></a>00277
- <a name="l00278"></a><a class="code" href="mpimod__multi_8f90.html#a4aceba15459fefd864a0ed3313b0073d">00278</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a4aceba15459fefd864a0ed3313b0073d">mpreadsp</a>(ktape,p,kdim,klev)
- <a name="l00279"></a>00279 <span class="keywordtype">real</span> p(kdim,klev)
- <a name="l00280"></a>00280 <span class="keyword">read</span> (ktape) p
- <a name="l00281"></a>00281 return
- <a name="l00282"></a>00282 <span class="keyword"> end</span>
- <a name="l00283"></a>00283
- <a name="l00284"></a><a class="code" href="mpimod__multi_8f90.html#a463456bde27045e2cf286e6e6082b9aa">00284</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a463456bde27045e2cf286e6e6082b9aa">mpreadgp</a>(ktape,p,kdim,klev)
- <a name="l00285"></a>00285 <span class="keywordtype">real</span> p(kdim,klev)
- <a name="l00286"></a>00286 <span class="keyword">read</span> (ktape) p
- <a name="l00287"></a>00287 return
- <a name="l00288"></a>00288 <span class="keyword"> end</span>
- <a name="l00289"></a>00289
- <a name="l00290"></a><a class="code" href="mpimod__multi_8f90.html#aca5ad2279542f783c1d862333da96744">00290</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#aca5ad2279542f783c1d862333da96744">mpwritesp</a>(ktape,p,kdim,klev)
- <a name="l00291"></a>00291 <span class="keywordtype">real</span> p(kdim,klev)
- <a name="l00292"></a>00292 <span class="keyword">write</span> (ktape) p
- <a name="l00293"></a>00293 return
- <a name="l00294"></a>00294 <span class="keyword"> end</span>
- <a name="l00295"></a>00295
- <a name="l00296"></a><a class="code" href="mpimod__multi_8f90.html#a3e3ab4b6cd8d7863f7bfc0e74b370488">00296</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a3e3ab4b6cd8d7863f7bfc0e74b370488">mpwritegp</a>(ktape,p,kdim,klev)
- <a name="l00297"></a>00297 <span class="keywordtype">real</span> p(kdim,klev)
- <a name="l00298"></a>00298 <span class="keyword">write</span> (ktape) p
- <a name="l00299"></a>00299 return
- <a name="l00300"></a>00300 <span class="keyword"> end</span>
- <a name="l00301"></a>00301
- <a name="l00302"></a><a class="code" href="mpimod__multi_8f90.html#a325e1b8f8412b422a06fb7558f212f7e">00302</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a325e1b8f8412b422a06fb7558f212f7e">mpwritegph</a>(ktape,p,kdim,klev,ihead)
- <a name="l00303"></a>00303 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00304"></a>00304 <span class="keywordtype">integer</span> :: ihead(8)
- <a name="l00305"></a>00305 <span class="keyword">write</span> (ktape) ihead
- <a name="l00306"></a>00306 <span class="keyword">write</span> (ktape) p
- <a name="l00307"></a>00307
- <a name="l00308"></a>00308 return
- <a name="l00309"></a>00309 <span class="keyword"> end</span>
- <a name="l00310"></a>00310
- <a name="l00311"></a>00311
- <a name="l00312"></a><a class="code" href="mpimod__multi_8f90.html#acb4faf87d9aa8c0bfc86d75f261989c5">00312</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#aa6e5f06de5db1764ec5ee24c21bae378">mpi_info</a>(nprocess,pid) ! get nproc and pid
- <a name="l00313"></a>00313 <span class="keywordtype">integer</span> nprocess, pid
- <a name="l00314"></a>00314 nprocess = 1
- <a name="l00315"></a>00315 pid = 0
- <a name="l00316"></a>00316 return
- <a name="l00317"></a>00317 <span class="keyword"> end subroutine mpi_info</span>
- <a name="l00318"></a>00318
- <a name="l00319"></a>00319
- <a name="l00320"></a><a class="code" href="mpimod__multi_8f90.html#acf82ae878fff75151cab59cdd0925ae0">00320</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#acf82ae878fff75151cab59cdd0925ae0">mpgetsp</a>(yn,p,kdim,klev)
- <a name="l00321"></a>00321 <span class="keywordtype">character (len=*)</span> :: yn
- <a name="l00322"></a>00322 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00323"></a>00323 call <a class="code" href="restartmod_8f90.html#af0f1ce9b6762aa2537cc22d5fc319b7c">get_restart_array</a>(yn,p,kdim,kdim,klev)
- <a name="l00324"></a>00324 return
- <a name="l00325"></a>00325 <span class="keyword"> end subroutine mpgetsp</span>
- <a name="l00326"></a>00326
- <a name="l00327"></a>00327
- <a name="l00328"></a><a class="code" href="mpimod__multi_8f90.html#a58d54c2e0590e63a7459417831afe5cf">00328</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a58d54c2e0590e63a7459417831afe5cf">mpgetgp</a>(yn,p,kdim,klev)
- <a name="l00329"></a>00329 <span class="keywordtype">character (len=*)</span> :: yn
- <a name="l00330"></a>00330 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00331"></a>00331 call <a class="code" href="restartmod_8f90.html#af0f1ce9b6762aa2537cc22d5fc319b7c">get_restart_array</a>(yn,p,kdim,kdim,klev)
- <a name="l00332"></a>00332 return
- <a name="l00333"></a>00333 <span class="keyword"> end subroutine mpgetgp</span>
- <a name="l00334"></a>00334
- <a name="l00335"></a>00335
- <a name="l00336"></a><a class="code" href="mpimod__multi_8f90.html#a79c341b7b52bf44470898581072660b8">00336</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a79c341b7b52bf44470898581072660b8">mpputsp</a>(yn,p,kdim,klev)
- <a name="l00337"></a>00337 <span class="keywordtype">character (len=*)</span> :: yn
- <a name="l00338"></a>00338 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00339"></a>00339 call <a class="code" href="restartmod_8f90.html#a52485001dbbaed032e48d894e6302c22">put_restart_array</a>(yn,p,kdim,kdim,klev)
- <a name="l00340"></a>00340 return
- <a name="l00341"></a>00341 <span class="keyword"> end subroutine mpputsp</span>
- <a name="l00342"></a>00342
- <a name="l00343"></a>00343
- <a name="l00344"></a><a class="code" href="mpimod__multi_8f90.html#a7e675330db7b46cf0bf0cc8edd2d413c">00344</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a7e675330db7b46cf0bf0cc8edd2d413c">mpputgp</a>(yn,p,kdim,klev)
- <a name="l00345"></a>00345 <span class="keywordtype">character (len=*)</span> :: yn
- <a name="l00346"></a>00346 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00347"></a>00347 call <a class="code" href="restartmod_8f90.html#a52485001dbbaed032e48d894e6302c22">put_restart_array</a>(yn,p,kdim,kdim,klev)
- <a name="l00348"></a>00348 return
- <a name="l00349"></a>00349 <span class="keyword"> end subroutine mpputgp</span>
- <a name="l00350"></a>00350
- <a name="l00351"></a>00351
- <a name="l00352"></a>00352 <span class="comment">! subroutine mpsurfgp(yn,p,kdim,klev)</span>
- <a name="l00353"></a>00353 <span class="comment">! character (len=*) :: yn</span>
- <a name="l00354"></a>00354 <span class="comment">! real :: p(kdim,klev)</span>
- <a name="l00355"></a>00355 <span class="comment">! call get_surf_array(yn,p,kdim,kdim,klev,iread)</span>
- <a name="l00356"></a>00356 <span class="comment">! return</span>
- <a name="l00357"></a>00357 <span class="comment">! end subroutine mpsurfgp</span>
- <a name="l00358"></a>00358 <span class="comment">!</span>
- <a name="l00359"></a>00359 <span class="comment">!</span>
- <a name="l00360"></a>00360 <span class="comment">! subroutine mpsurfyear(yn,p,kdim,kmon)</span>
- <a name="l00361"></a>00361 <span class="comment">! character (len=*) :: yn</span>
- <a name="l00362"></a>00362 <span class="comment">! real :: p(kdim,kmon)</span>
- <a name="l00363"></a>00363 <span class="comment">! call get_surf_year(yn,p,kdim,kmon,iread)</span>
- <a name="l00364"></a>00364 <span class="comment">! return</span>
- <a name="l00365"></a>00365 <span class="comment">! end subroutine mpsurfyear</span>
- <a name="l00366"></a>00366 <span class="comment">!</span>
- <a name="l00367"></a>00367 <span class="comment">!</span>
- <a name="l00368"></a>00368 <span class="comment">! subroutine mp3dyear(yn,p,kdim,klev,kmon)</span>
- <a name="l00369"></a>00369 <span class="comment">! character (len=*) :: yn</span>
- <a name="l00370"></a>00370 <span class="comment">! real :: p(kdim,klev,kmon)</span>
- <a name="l00371"></a>00371 <span class="comment">! call get_3d_year(yn,p,kdim,klev,kmon,iread)</span>
- <a name="l00372"></a>00372 <span class="comment">! return</span>
- <a name="l00373"></a>00373 <span class="comment">! end subroutine mp3dyear</span>
- <a name="l00374"></a>00374
- <a name="l00375"></a>00375
- <a name="l00376"></a><a class="code" href="mpimod__multi_8f90.html#a1b6ac2b98059a43359ac0edfeb9c2ad7">00376</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#a1b6ac2b98059a43359ac0edfeb9c2ad7">mpmaxval</a>(p,kdim,klev,pmax)
- <a name="l00377"></a>00377 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00378"></a>00378 pmax = maxval(p(:,:))
- <a name="l00379"></a>00379 return
- <a name="l00380"></a>00380 <span class="keyword"> end subroutine mpmaxval</span>
- <a name="l00381"></a>00381
- <a name="l00382"></a>00382
- <a name="l00383"></a><a class="code" href="mpimod__multi_8f90.html#ac1dfb34daad89cf72ff04b6a58919b2b">00383</a> <span class="keyword">subroutine </span><a class="code" href="mpimod_8f90.html#ac1dfb34daad89cf72ff04b6a58919b2b">mpsumval</a>(p,kdim,klev,psum)
- <a name="l00384"></a>00384 <span class="keywordtype">real</span> :: p(kdim,klev)
- <a name="l00385"></a>00385 psum = sum(p(:,:))
- <a name="l00386"></a>00386 return
- <a name="l00387"></a>00387 <span class="keyword"> end subroutine mpsumval</span>
- <a name="l00388"></a>00388
- </pre></div></div>
- </div>
- <div id="nav-path" class="navpath">
- <ul>
- <li class="navelem"><a class="el" href="mpimod__multi_8f90.html">mpimod_multi.f90</a> </li>
- <!-- window showing the filter options -->
- <div id="MSearchSelectWindow"
- onmouseover="return searchBox.OnSearchSelectShow()"
- onmouseout="return searchBox.OnSearchSelectHide()"
- onkeydown="return searchBox.OnSearchSelectKey(event)">
- <a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Defines</a></div>
- <!-- iframe showing the search results (closed by default) -->
- <div id="MSearchResultsWindow">
- <iframe src="javascript:void(0)" frameborder="0"
- name="MSearchResults" id="MSearchResults">
- </iframe>
- </div>
- <li class="footer">Generated on Wed Sep 21 2011 12:35:46 for PUMA by
- <a href="http://www.doxygen.org/index.html">
- <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.5.1 </li>
- </ul>
- </div>
- </body>
- </html>
|