qmpi.F90 57 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030203120322033203420352036203720382039204020412042204320442045204620472048204920502051205220532054205520562057205820592060206120622063206420652066206720682069207020712072
  1. #if defined(QMPI)
  2. module qmpi
  3. !
  4. ! A module defining a minimalist interface to a subset of MPI.
  5. ! The first five primitives can in theory be used to parallelize
  6. ! any program. The module hides type specification, communicators,
  7. ! explicit error handling, the need to give explicit buffer size etc.
  8. ! Also provided are a few interfaces for often used broadcast and
  9. ! reduction operations
  10. !
  11. ! © Helge Avlesen <avle@ii.uib.no>, para//ab
  12. !
  13. ! primitives: (optional arguments in brackets)
  14. !
  15. ! subroutine start_mpi()
  16. ! starts the mpi subsystem. all processesors are assigned a number (myid).
  17. ! the number of processors is numproc.
  18. ! subroutine stop_mpi()
  19. ! stops the mpi subsystem
  20. ! subroutine barrier([label])
  21. ! syncronization point for all processors. optionally prints a label on
  22. ! the master processor (0).
  23. ! subroutine send(data, target [,tag])
  24. ! send object data to processor number target, tag is an optional integer
  25. ! that defaults to 0. (if multiple messages are exchanged between a
  26. ! pair of processors, a unique tag must be used for each exhange)
  27. ! subroutine receive(data, source [,tag])
  28. ! get object data from processor source, tag is optional and as for send
  29. ! MPI will fail if the size of the object received is different from what
  30. ! was sent.
  31. !
  32. ! The rest of the routines are included for convenience, they can be
  33. ! also be implemented using the above subroutines.
  34. !
  35. ! subroutine broadcast(data [,root])
  36. ! broadcast data (any type) from processor root (default=0) to all
  37. ! other processors.
  38. ! subroutine mbroadcast(data [,data2,data3,data4,data5,data6] [,root])
  39. ! broadcast up to 6 scalar variables of the same type, to all processors
  40. ! from processor root (default=0)
  41. ! subroutine reduce(type, data [,data2,data3,data4,data5,data6] [,root] )
  42. ! reduce the scalar data, optionally also data2-data6, return result
  43. ! on all processes. the operation can currently be of type 'sum', 'mul',
  44. ! 'min' or 'max' i.e. a sum or a product. data-data6 must be of the
  45. ! same type. if integer root is present, only return result on that
  46. ! processor (faster)
  47. !
  48. ! Example: a program that sends a real from processor 0 to processor 1
  49. ! use qmpi
  50. ! real data
  51. ! call start_mpi
  52. ! data=myid
  53. ! if(myid==0) call send(data, 1)
  54. ! if(myid==1) then
  55. ! call receive(data, 0)
  56. ! print *,'hello, I am',myid,'got ',data,'from process 0'
  57. ! end if
  58. ! call stop_mpi
  59. ! end
  60. !
  61. ! More advanced usage example: to send a derived type from 0 to 1;
  62. ! pack it in a string (could be packed into any array), send, receive, unpack.
  63. !
  64. ! type(any_type) var1
  65. ! character, allocatable :: buffer(:)
  66. ! ...
  67. ! N=size(transfer(var1,(/'x'/)))) !! compute size of type once
  68. ! allocate(buffer(N))
  69. ! if(myid==0)then
  70. ! buffer = transfer(var1,buffer)
  71. ! call send(buffer,1)
  72. ! end if
  73. ! if(myid==1)then
  74. ! call receive(buffer,0)
  75. ! var1 = transfer(buffer,var1)
  76. ! end if
  77. ! ...
  78. !
  79. #warning "COMPILING WITH QMPI CODE"
  80. include 'mpif.h'
  81. integer, public :: qmpi_proc_num, qmpi_num_proc, ierr, errorcode, mpistatus(mpi_status_size)
  82. logical, public :: master=.false., slave=.false.
  83. ! some kinds. could use selected_real_kind(..) for this instead of hard coding
  84. integer, parameter :: dp=8, sp=4, long=8, short=2
  85. interface send
  86. module procedure &
  87. qmpi_send_real4, &
  88. qmpi_send_real4_1d, &
  89. qmpi_send_real4_2d, &
  90. qmpi_send_real4_3d, &
  91. qmpi_send_real4_4d, &
  92. qmpi_send_real8, &
  93. qmpi_send_real8_1d, &
  94. qmpi_send_real8_2d, &
  95. qmpi_send_real8_3d, &
  96. qmpi_send_real8_4d, &
  97. qmpi_send_integer4, &
  98. qmpi_send_integer4_1d, &
  99. qmpi_send_integer4_2d, &
  100. qmpi_send_integer4_3d, &
  101. qmpi_send_integer4_4d, &
  102. qmpi_send_integer8, &
  103. qmpi_send_integer8_1d, &
  104. qmpi_send_integer8_2d, &
  105. qmpi_send_integer8_3d, &
  106. qmpi_send_integer8_4d, &
  107. qmpi_send_string, &
  108. qmpi_send_character_1d,&
  109. qmpi_send_logical
  110. end interface
  111. interface receive
  112. module procedure &
  113. qmpi_recv_real4, &
  114. qmpi_recv_real4_1d, &
  115. qmpi_recv_real4_2d, &
  116. qmpi_recv_real4_3d, &
  117. qmpi_recv_real4_4d, &
  118. qmpi_recv_real8, &
  119. qmpi_recv_real8_1d, &
  120. qmpi_recv_real8_2d, &
  121. qmpi_recv_real8_3d, &
  122. qmpi_recv_real8_4d, &
  123. qmpi_recv_integer4, &
  124. qmpi_recv_integer4_1d, &
  125. qmpi_recv_integer4_2d, &
  126. qmpi_recv_integer4_3d, &
  127. qmpi_recv_integer4_4d, &
  128. qmpi_recv_integer8, &
  129. qmpi_recv_integer8_1d, &
  130. qmpi_recv_integer8_2d, &
  131. qmpi_recv_integer8_3d, &
  132. qmpi_recv_integer8_4d, &
  133. qmpi_recv_string, &
  134. qmpi_recv_character_1d,&
  135. qmpi_recv_logical
  136. end interface
  137. interface reduce
  138. module procedure &
  139. qmpi_integer_reduction, &
  140. qmpi_integer8_reduction,&
  141. qmpi_real_reduction, &
  142. qmpi_real8_reduction
  143. end interface
  144. interface broadcast
  145. module procedure &
  146. qmpi_broadcast_logical, &
  147. qmpi_broadcast_string, &
  148. qmpi_broadcast_stringarr,&
  149. qmpi_broadcast_integer4, &
  150. qmpi_broadcast_integer4_array1d, &
  151. qmpi_broadcast_integer4_array2d, &
  152. qmpi_broadcast_integer8, &
  153. qmpi_broadcast_integer8_array1d, &
  154. qmpi_broadcast_integer8_array2d, &
  155. qmpi_broadcast_real4, &
  156. qmpi_broadcast_real4_array1d, &
  157. qmpi_broadcast_real4_array2d, &
  158. qmpi_broadcast_real4_array3d, &
  159. qmpi_broadcast_real4_array4d, &
  160. qmpi_broadcast_real8, &
  161. qmpi_broadcast_real8_array1d, &
  162. qmpi_broadcast_real8_array2d, &
  163. qmpi_broadcast_real8_array3d, &
  164. qmpi_broadcast_real8_array4d
  165. end interface
  166. interface mbroadcast
  167. module procedure &
  168. qmpi_broadcast_logicals, &
  169. qmpi_broadcast_real4s, &
  170. qmpi_broadcast_real8s, &
  171. qmpi_broadcast_integer4s, &
  172. qmpi_broadcast_integer8s
  173. end interface
  174. contains
  175. subroutine start_mpi()
  176. !
  177. ! initialize the core MPI subsystem
  178. ! this routine should be called as the first statement in the program.
  179. ! MPI does not specify what happen before MPI_init and after mpi_finalize
  180. !
  181. implicit none
  182. call mpi_init(ierr)
  183. call mpi_comm_size(mpi_comm_world, qmpi_num_proc, ierr)
  184. call mpi_comm_rank(mpi_comm_world, qmpi_proc_num, ierr)
  185. master=.false.
  186. if(qmpi_proc_num==0) master=.true.
  187. if(qmpi_proc_num>0) slave=.true.
  188. print*,'Inne i start_mpi: qmpi_proc_num =',qmpi_proc_num,' master =',master
  189. if(master) then
  190. write(*,'(a,i0,a)') 'MPI started with ',qmpi_num_proc,' processors'
  191. end if
  192. end subroutine start_mpi
  193. subroutine stop_mpi()
  194. implicit none
  195. call mpi_finalize(ierr)
  196. stop
  197. end subroutine stop_mpi
  198. subroutine barrier(label)
  199. ! makes all processes sync at this point, optionally print a label
  200. implicit none
  201. character(*), optional :: label
  202. call mpi_barrier(mpi_comm_world, ierr)
  203. if(master.and.present(label)) print *,'---barrier---',label,'---------'
  204. end subroutine barrier
  205. subroutine qmpi_send_logical(data, target, tag)
  206. implicit none
  207. logical data
  208. integer target, counter, given_tag
  209. integer, optional :: tag
  210. given_tag=0
  211. if(present(tag)) given_tag=tag
  212. counter=1
  213. call mpi_send(data, counter, mpi_logical, target, given_tag, mpi_comm_world, ierr)
  214. if(ierr.ne.0)then
  215. print *,'error send_logical count=',counter,'tag=',given_tag
  216. stop
  217. end if
  218. end subroutine qmpi_send_logical
  219. subroutine qmpi_send_string(data, target, tag)
  220. implicit none
  221. character(*) data
  222. integer target, counter, given_tag
  223. integer, optional :: tag
  224. given_tag=0
  225. if(present(tag)) given_tag=tag
  226. counter=len(data)
  227. call mpi_send(data, counter, mpi_character, target, given_tag, mpi_comm_world, ierr)
  228. if(ierr.ne.0)then
  229. print *,'error send_string count=',counter,'tag=',given_tag
  230. stop
  231. end if
  232. end subroutine qmpi_send_string
  233. subroutine qmpi_send_character_1d(data, target, tag)
  234. implicit none
  235. character data(:)
  236. integer target, counter, given_tag
  237. integer, optional :: tag
  238. given_tag=0
  239. if(present(tag)) given_tag=tag
  240. counter=size(data)
  241. call mpi_send(data, counter, mpi_character, target, given_tag, mpi_comm_world, ierr)
  242. if(ierr.ne.0)then
  243. print *,'error send_character_1d count=',counter,'tag=',given_tag
  244. stop
  245. end if
  246. end subroutine qmpi_send_character_1d
  247. subroutine qmpi_recv_character_1d(data, target, tag)
  248. implicit none
  249. character data(:)
  250. integer target, counter, given_tag
  251. integer, optional :: tag
  252. given_tag=0
  253. if(present(tag)) given_tag=tag
  254. counter=size(data)
  255. call mpi_recv(data, counter, mpi_character, target, given_tag, mpi_comm_world, mpistatus, ierr)
  256. if(ierr.ne.0)then
  257. print *,'error recv_character_1d count=',counter,'tag=',given_tag
  258. stop
  259. end if
  260. end subroutine qmpi_recv_character_1d
  261. subroutine qmpi_send_integer4(data, target, tag)
  262. implicit none
  263. integer(sp) data
  264. integer target, counter, given_tag
  265. integer, optional :: tag
  266. given_tag=0
  267. if(present(tag)) given_tag=tag
  268. counter=1
  269. call mpi_send(data, counter, mpi_integer, target, given_tag, mpi_comm_world, ierr)
  270. if(ierr.ne.0)then
  271. print *,'error send_integer4 count=',counter,'tag=',given_tag
  272. stop
  273. end if
  274. end subroutine qmpi_send_integer4
  275. subroutine qmpi_send_integer4_1d(data, target, tag)
  276. implicit none
  277. integer(sp) data(:)
  278. integer target, counter, given_tag
  279. integer, optional :: tag
  280. given_tag=0
  281. if(present(tag)) given_tag=tag
  282. counter=size(data)
  283. call mpi_send(data, counter, mpi_integer, target, given_tag, mpi_comm_world, ierr)
  284. if(ierr.ne.0)then
  285. print *,'error send_integer4_1d count=',counter,'tag=',given_tag
  286. stop
  287. end if
  288. end subroutine qmpi_send_integer4_1d
  289. subroutine qmpi_send_integer4_2d(data, target, tag)
  290. implicit none
  291. integer(sp) data(:,:)
  292. integer target, counter, given_tag
  293. integer, optional :: tag
  294. given_tag=0
  295. if(present(tag)) given_tag=tag
  296. counter=size(data,1)*size(data,2)
  297. call mpi_send(data, counter, mpi_integer, target, given_tag, mpi_comm_world, ierr)
  298. if(ierr.ne.0)then
  299. print *,'error send_integer4_2d count=',counter,'tag=',given_tag
  300. stop
  301. end if
  302. end subroutine qmpi_send_integer4_2d
  303. subroutine qmpi_send_integer4_3d(data, target, tag)
  304. implicit none
  305. integer(sp) data(:,:,:)
  306. integer target, counter, given_tag
  307. integer, optional :: tag
  308. given_tag=0
  309. if(present(tag)) given_tag=tag
  310. counter=size(data,1)*size(data,2)*size(data,3)
  311. call mpi_send(data, counter, mpi_integer, target, given_tag, mpi_comm_world, ierr)
  312. if(ierr.ne.0)then
  313. print *,'error send_integer4_3d count=',counter,'tag=',given_tag
  314. stop
  315. end if
  316. end subroutine qmpi_send_integer4_3d
  317. subroutine qmpi_send_integer4_4d(data, target, tag)
  318. implicit none
  319. integer(sp) data(:,:,:,:)
  320. integer target, counter, given_tag
  321. integer, optional :: tag
  322. given_tag=0
  323. if(present(tag)) given_tag=tag
  324. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  325. call mpi_send(data, counter, mpi_integer, target, given_tag, mpi_comm_world, ierr)
  326. if(ierr.ne.0)then
  327. print *,'error send_integer4_4d count=',counter,'tag=',given_tag
  328. stop
  329. end if
  330. end subroutine qmpi_send_integer4_4d
  331. subroutine qmpi_send_integer8(data, target, tag)
  332. implicit none
  333. integer(long) data
  334. integer target, counter, given_tag
  335. integer, optional :: tag
  336. given_tag=0
  337. if(present(tag)) given_tag=tag
  338. counter=1
  339. call mpi_send(data, counter, mpi_integer8, target, given_tag, mpi_comm_world, ierr)
  340. if(ierr.ne.0)then
  341. print *,'error send_integer8 count=',counter,'tag=',given_tag
  342. stop
  343. end if
  344. end subroutine qmpi_send_integer8
  345. subroutine qmpi_send_integer8_1d(data, target, tag)
  346. implicit none
  347. integer(long) data(:)
  348. integer target, counter, given_tag
  349. integer, optional :: tag
  350. given_tag=0
  351. if(present(tag)) given_tag=tag
  352. counter=size(data)
  353. call mpi_send(data, counter, mpi_integer8, target, given_tag, mpi_comm_world, ierr)
  354. if(ierr.ne.0)then
  355. print *,'error send_integer8_1d count=',counter,'tag=',given_tag
  356. stop
  357. end if
  358. end subroutine qmpi_send_integer8_1d
  359. subroutine qmpi_send_integer8_2d(data, target, tag)
  360. implicit none
  361. integer(long) data(:,:)
  362. integer target, counter, given_tag
  363. integer, optional :: tag
  364. given_tag=0
  365. if(present(tag)) given_tag=tag
  366. counter=size(data,1)*size(data,2)
  367. call mpi_send(data, counter, mpi_integer8, target, given_tag, mpi_comm_world, ierr)
  368. if(ierr.ne.0)then
  369. print *,'error send_integer8_2d count=',counter,'tag=',given_tag
  370. stop
  371. end if
  372. end subroutine qmpi_send_integer8_2d
  373. subroutine qmpi_send_integer8_3d(data, target, tag)
  374. implicit none
  375. integer(8) data(:,:,:)
  376. integer target, counter, given_tag
  377. integer, optional :: tag
  378. given_tag=0
  379. if(present(tag)) given_tag=tag
  380. counter=size(data,1)*size(data,2)*size(data,3)
  381. call mpi_send(data, counter, mpi_integer8, target, given_tag, mpi_comm_world, ierr)
  382. if(ierr.ne.0)then
  383. print *,'error send_integer8_3d count=',counter,'tag=',given_tag
  384. stop
  385. end if
  386. end subroutine qmpi_send_integer8_3d
  387. subroutine qmpi_send_integer8_4d(data, target, tag)
  388. implicit none
  389. integer(8) data(:,:,:,:)
  390. integer target, counter, given_tag
  391. integer, optional :: tag
  392. given_tag=0
  393. if(present(tag)) given_tag=tag
  394. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  395. call mpi_send(data, counter, mpi_integer8, target, given_tag, mpi_comm_world, ierr)
  396. if(ierr.ne.0)then
  397. print *,'error send_integer8_4d count=',counter,'tag=',given_tag
  398. stop
  399. end if
  400. end subroutine qmpi_send_integer8_4d
  401. subroutine qmpi_send_real4(data, target, tag)
  402. implicit none
  403. real(sp) data
  404. integer target
  405. integer, optional :: tag
  406. integer counter, given_tag
  407. given_tag=0
  408. if(present(tag)) given_tag=tag
  409. counter=1
  410. call mpi_send(data, counter, mpi_real, target, given_tag, mpi_comm_world, ierr)
  411. if(ierr.ne.0)then
  412. print *,'error send_real4 count=',counter,'tag=',given_tag
  413. stop
  414. end if
  415. end subroutine qmpi_send_real4
  416. subroutine qmpi_send_real8(data, target, tag)
  417. implicit none
  418. real(dp) data
  419. integer target
  420. integer, optional :: tag
  421. integer counter, given_tag
  422. given_tag=0
  423. if(present(tag)) given_tag=tag
  424. counter=1
  425. call mpi_send(data, counter, mpi_double_precision, target, given_tag, mpi_comm_world, ierr)
  426. if(ierr.ne.0)then
  427. print *,'error send_real8 count=',counter,'tag=',given_tag
  428. stop
  429. end if
  430. end subroutine qmpi_send_real8
  431. subroutine qmpi_send_real4_1d(data, target, tag)
  432. implicit none
  433. real(sp) data(:)
  434. integer target
  435. integer, optional :: tag
  436. integer counter, given_tag
  437. given_tag=0
  438. if(present(tag)) given_tag=tag
  439. counter=size(data)
  440. call mpi_send(data, counter, mpi_real, target, given_tag, mpi_comm_world, ierr)
  441. if(ierr.ne.0)then
  442. print *,'error send_real4_1d count=',counter,'tag=',given_tag
  443. stop
  444. end if
  445. end subroutine qmpi_send_real4_1d
  446. subroutine qmpi_send_real8_1d(data, target, tag)
  447. implicit none
  448. real(dp) data(:)
  449. integer target
  450. integer, optional :: tag
  451. integer counter, given_tag
  452. given_tag=0
  453. if(present(tag)) given_tag=tag
  454. counter=size(data)
  455. call mpi_send(data, counter, mpi_double_precision, target, given_tag, mpi_comm_world, ierr)
  456. if(ierr.ne.0)then
  457. print *,'error send_real8_1d count=',counter,'tag=',given_tag
  458. stop
  459. end if
  460. end subroutine qmpi_send_real8_1d
  461. subroutine qmpi_send_real4_2d(data, target, tag)
  462. implicit none
  463. real(sp) data(:,:)
  464. integer target
  465. integer, optional :: tag
  466. integer counter, given_tag
  467. given_tag=0
  468. if(present(tag)) given_tag=tag
  469. counter=size(data,1)*size(data,2)
  470. call mpi_send(data, counter, mpi_real, target, given_tag, mpi_comm_world, ierr)
  471. if(ierr.ne.0)then
  472. print *,'error send_real4_2d count=',counter,'tag=',given_tag
  473. stop
  474. end if
  475. end subroutine qmpi_send_real4_2d
  476. subroutine qmpi_send_real8_2d(data, target, tag)
  477. implicit none
  478. real(dp) data(:,:)
  479. integer target
  480. integer, optional :: tag
  481. integer counter, given_tag
  482. given_tag=0
  483. if(present(tag)) given_tag=tag
  484. counter=size(data,1)*size(data,2)
  485. call mpi_send(data, counter, mpi_double_precision, target, given_tag, mpi_comm_world, ierr)
  486. if(ierr.ne.0)then
  487. print *,'error send_real8_2d count=',counter,'tag=',given_tag
  488. stop
  489. end if
  490. end subroutine qmpi_send_real8_2d
  491. subroutine qmpi_send_real4_3d(data, target, tag)
  492. implicit none
  493. real(sp) data(:,:,:)
  494. integer target
  495. integer, optional :: tag
  496. integer counter, given_tag
  497. given_tag=0
  498. if(present(tag)) given_tag=tag
  499. counter=size(data,1)*size(data,2)*size(data,3)
  500. call mpi_send(data, counter, mpi_real, target, given_tag, mpi_comm_world, ierr)
  501. if(ierr.ne.0)then
  502. print *,'error send_real4_3d count=',counter,'tag=',given_tag
  503. stop
  504. end if
  505. end subroutine qmpi_send_real4_3d
  506. subroutine qmpi_send_real8_3d(data, target, tag)
  507. implicit none
  508. real(dp) data(:,:,:)
  509. integer target
  510. integer, optional :: tag
  511. integer counter, given_tag
  512. given_tag=0
  513. if(present(tag)) given_tag=tag
  514. counter=size(data,1)*size(data,2)*size(data,3)
  515. call mpi_send(data, counter, mpi_double_precision, target, given_tag, mpi_comm_world, ierr)
  516. if(ierr.ne.0)then
  517. print *,'error send_real8_3d count=',counter,'tag=',given_tag
  518. stop
  519. end if
  520. end subroutine qmpi_send_real8_3d
  521. subroutine qmpi_send_real4_4d(data, target, tag)
  522. implicit none
  523. real(sp) data(:,:,:,:)
  524. integer target
  525. integer, optional :: tag
  526. integer counter, given_tag
  527. given_tag=0
  528. if(present(tag)) given_tag=tag
  529. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  530. call mpi_send(data, counter, mpi_real, target, given_tag, mpi_comm_world, ierr)
  531. if(ierr.ne.0)then
  532. print *,'error send_real4_4d count=',counter,'tag=',given_tag
  533. stop
  534. end if
  535. end subroutine qmpi_send_real4_4d
  536. subroutine qmpi_send_real8_4d(data, target, tag)
  537. implicit none
  538. real(dp) data(:,:,:,:)
  539. integer target
  540. integer, optional :: tag
  541. integer counter, given_tag
  542. given_tag=0
  543. if(present(tag)) given_tag=tag
  544. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  545. call mpi_send(data, counter, mpi_double_precision, target, given_tag, mpi_comm_world, ierr)
  546. if(ierr.ne.0)then
  547. print *,'error send_real8_4d count=',counter,'tag=',given_tag
  548. stop
  549. end if
  550. end subroutine qmpi_send_real8_4d
  551. subroutine qmpi_recv_integer4(data, source, tag)
  552. implicit none
  553. integer(sp) data
  554. integer source, counter, given_tag
  555. integer, optional :: tag
  556. given_tag=0
  557. if(present(tag)) given_tag=tag
  558. counter=1
  559. call mpi_recv(data, counter, mpi_integer, source, given_tag, mpi_comm_world, mpistatus, ierr)
  560. if(ierr.ne.0)then
  561. print *,'error recv_integer4_1d count=',counter,'tag=',given_tag
  562. stop
  563. end if
  564. end subroutine qmpi_recv_integer4
  565. subroutine qmpi_recv_integer4_1d(data, source, tag)
  566. implicit none
  567. integer(sp) data(:)
  568. integer source, counter, given_tag
  569. integer, optional :: tag
  570. given_tag=0
  571. if(present(tag)) given_tag=tag
  572. counter=size(data)
  573. call mpi_recv(data, counter, mpi_integer, source, given_tag, mpi_comm_world, mpistatus, ierr)
  574. if(ierr.ne.0)then
  575. print *,'error recv_integer4_1d count=',counter,'tag=',given_tag
  576. stop
  577. end if
  578. end subroutine qmpi_recv_integer4_1d
  579. subroutine qmpi_recv_integer4_2d(data, source, tag)
  580. implicit none
  581. integer(sp) data(:,:)
  582. integer source, counter, given_tag
  583. integer, optional :: tag
  584. given_tag=0
  585. if(present(tag)) given_tag=tag
  586. counter=size(data,1)*size(data,2)
  587. call mpi_recv(data, counter, mpi_integer, source, given_tag, mpi_comm_world, mpistatus, ierr)
  588. if(ierr.ne.0)then
  589. print *,'error recv_integer4_2d count=',counter,'tag=',given_tag
  590. stop
  591. end if
  592. end subroutine qmpi_recv_integer4_2d
  593. subroutine qmpi_recv_integer4_3d(data, source, tag)
  594. implicit none
  595. integer(sp) data(:,:,:)
  596. integer source, counter, given_tag
  597. integer, optional :: tag
  598. given_tag=0
  599. if(present(tag)) given_tag=tag
  600. counter=size(data,1)*size(data,2)*size(data,3)
  601. call mpi_recv(data, counter, mpi_integer, source, given_tag, mpi_comm_world, mpistatus, ierr)
  602. if(ierr.ne.0)then
  603. print *,'error recv_integer4_3d count=',counter,'tag=',given_tag
  604. stop
  605. end if
  606. end subroutine qmpi_recv_integer4_3d
  607. subroutine qmpi_recv_integer4_4d(data, source, tag)
  608. implicit none
  609. integer(sp) data(:,:,:,:)
  610. integer source, counter, given_tag
  611. integer, optional :: tag
  612. given_tag=0
  613. if(present(tag)) given_tag=tag
  614. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  615. call mpi_recv(data, counter, mpi_integer, source, given_tag, mpi_comm_world, mpistatus, ierr)
  616. if(ierr.ne.0)then
  617. print *,'error recv_integer4_4d count=',counter,'tag=',given_tag
  618. stop
  619. end if
  620. end subroutine qmpi_recv_integer4_4d
  621. subroutine qmpi_recv_integer8(data, source, tag)
  622. implicit none
  623. integer(long) data
  624. integer source, counter, given_tag
  625. integer, optional :: tag
  626. given_tag=0
  627. if(present(tag)) given_tag=tag
  628. counter=1
  629. call mpi_recv(data, counter, mpi_integer8, source, given_tag, mpi_comm_world, mpistatus, ierr)
  630. if(ierr.ne.0)then
  631. print *,'error recv_integer8 count=',counter,'tag=',given_tag
  632. stop
  633. end if
  634. end subroutine qmpi_recv_integer8
  635. subroutine qmpi_recv_integer8_1d(data, source, tag)
  636. implicit none
  637. integer(long) data(:)
  638. integer source, counter, given_tag
  639. integer, optional :: tag
  640. given_tag=0
  641. if(present(tag)) given_tag=tag
  642. counter=size(data)
  643. call mpi_recv(data, counter, mpi_integer8, source, given_tag, mpi_comm_world, mpistatus, ierr)
  644. if(ierr.ne.0)then
  645. print *,'error recv_integer8_1d count=',counter,'tag=',given_tag
  646. stop
  647. end if
  648. end subroutine qmpi_recv_integer8_1d
  649. subroutine qmpi_recv_integer8_2d(data, source, tag)
  650. implicit none
  651. integer(long) data(:,:)
  652. integer source, counter, given_tag
  653. integer, optional :: tag
  654. given_tag=0
  655. if(present(tag)) given_tag=tag
  656. counter=size(data,1)*size(data,2)
  657. call mpi_recv(data, counter, mpi_integer8, source, given_tag, mpi_comm_world, mpistatus, ierr)
  658. if(ierr.ne.0)then
  659. print *,'error recv_integer8_2d count=',counter,'tag=',given_tag
  660. stop
  661. end if
  662. end subroutine qmpi_recv_integer8_2d
  663. subroutine qmpi_recv_integer8_3d(data, source, tag)
  664. implicit none
  665. integer(8) data(:,:,:)
  666. integer source, counter, given_tag
  667. integer, optional :: tag
  668. given_tag=0
  669. if(present(tag)) given_tag=tag
  670. counter=size(data,1)*size(data,2)*size(data,3)
  671. call mpi_recv(data, counter, mpi_integer8, source, given_tag, mpi_comm_world, mpistatus, ierr)
  672. if(ierr.ne.0)then
  673. print *,'error recv_integer8_3d count=',counter,'tag=',given_tag
  674. stop
  675. end if
  676. end subroutine qmpi_recv_integer8_3d
  677. subroutine qmpi_recv_integer8_4d(data, source, tag)
  678. implicit none
  679. integer(8) data(:,:,:,:)
  680. integer source, counter, given_tag
  681. integer, optional :: tag
  682. given_tag=0
  683. if(present(tag)) given_tag=tag
  684. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  685. call mpi_recv(data, counter, mpi_integer8, source, given_tag, mpi_comm_world, mpistatus, ierr)
  686. if(ierr.ne.0)then
  687. print *,'error recv_integer8_4d count=',counter,'tag=',given_tag
  688. stop
  689. end if
  690. end subroutine qmpi_recv_integer8_4d
  691. subroutine qmpi_recv_real4(data, source, tag)
  692. implicit none
  693. real(sp) data
  694. integer source
  695. integer, optional :: tag
  696. integer counter, given_tag
  697. given_tag=0
  698. if(present(tag)) given_tag=tag
  699. counter=1
  700. call mpi_recv(data, counter, mpi_real, source, given_tag, mpi_comm_world, mpistatus, ierr)
  701. if(ierr.ne.0)then
  702. print *,'error recv_real4 count=',counter,'tag=',given_tag
  703. stop
  704. end if
  705. end subroutine qmpi_recv_real4
  706. subroutine qmpi_recv_real8(data, source, tag)
  707. implicit none
  708. real(dp) data
  709. integer source
  710. integer, optional :: tag
  711. integer counter, given_tag
  712. given_tag=0
  713. if(present(tag)) given_tag=tag
  714. counter=1
  715. call mpi_recv(data, counter, mpi_double_precision, source, given_tag, mpi_comm_world, mpistatus, ierr)
  716. if(ierr.ne.0)then
  717. print *,'error recv_real8 count=',counter,'tag=',given_tag
  718. stop
  719. end if
  720. end subroutine qmpi_recv_real8
  721. subroutine qmpi_recv_real4_1d(data, source, tag)
  722. implicit none
  723. real(sp) data(:)
  724. integer source
  725. integer, optional :: tag
  726. integer counter, given_tag
  727. given_tag=0
  728. if(present(tag)) given_tag=tag
  729. counter=size(data)
  730. call mpi_recv(data, counter, mpi_real, source, given_tag, mpi_comm_world, mpistatus, ierr)
  731. if(ierr.ne.0)then
  732. print *,'error recv_real4_1d count=',counter,'tag=',given_tag
  733. stop
  734. end if
  735. end subroutine qmpi_recv_real4_1d
  736. subroutine qmpi_recv_real8_1d(data, source, tag)
  737. implicit none
  738. real(dp) data(:)
  739. integer source
  740. integer, optional :: tag
  741. integer counter, given_tag
  742. given_tag=0
  743. if(present(tag)) given_tag=tag
  744. counter=size(data)
  745. call mpi_recv(data, counter, mpi_double_precision, source, given_tag, mpi_comm_world, mpistatus, ierr)
  746. if(ierr.ne.0)then
  747. print *,'error recv_real8_1d count=',counter,'tag=',given_tag
  748. stop
  749. end if
  750. end subroutine qmpi_recv_real8_1d
  751. subroutine qmpi_recv_real4_2d(data, source, tag)
  752. implicit none
  753. real(sp) data(:,:)
  754. integer source
  755. integer, optional :: tag
  756. integer counter, given_tag
  757. given_tag=0
  758. if(present(tag)) given_tag=tag
  759. counter=size(data,1)*size(data,2)
  760. call mpi_recv(data, counter, mpi_real, source, given_tag, mpi_comm_world, mpistatus, ierr)
  761. if(ierr.ne.0)then
  762. print *,'error recv_real4_2d count=',counter,'tag=',given_tag
  763. stop
  764. end if
  765. end subroutine qmpi_recv_real4_2d
  766. subroutine qmpi_recv_real8_2d(data, source, tag)
  767. implicit none
  768. real(dp) data(:,:)
  769. integer source
  770. integer, optional :: tag
  771. integer counter, given_tag
  772. given_tag=0
  773. if(present(tag)) given_tag=tag
  774. counter=size(data,1)*size(data,2)
  775. call mpi_recv(data, counter, mpi_double_precision, source, given_tag, mpi_comm_world, mpistatus, ierr)
  776. if(ierr.ne.0)then
  777. print *,'error recv_real8_2d count=',counter,'tag=',given_tag
  778. stop
  779. end if
  780. end subroutine qmpi_recv_real8_2d
  781. subroutine qmpi_recv_real4_3d(data, source, tag)
  782. implicit none
  783. real(sp) data(:,:,:)
  784. integer source
  785. integer, optional :: tag
  786. integer counter, given_tag
  787. given_tag=0
  788. if(present(tag)) given_tag=tag
  789. counter=size(data,1)*size(data,2)*size(data,3)
  790. call mpi_recv(data, counter, mpi_real, source, given_tag, mpi_comm_world, mpistatus, ierr)
  791. if(ierr.ne.0)then
  792. print *,'error recv_real4_3d count=',counter,'tag=',given_tag
  793. stop
  794. end if
  795. end subroutine qmpi_recv_real4_3d
  796. subroutine qmpi_recv_real8_3d(data, source, tag)
  797. implicit none
  798. real(dp) data(:,:,:)
  799. integer source
  800. integer, optional :: tag
  801. integer counter, given_tag
  802. given_tag=0
  803. if(present(tag)) given_tag=tag
  804. counter=size(data,1)*size(data,2)*size(data,3)
  805. call mpi_recv(data, counter, mpi_double_precision, source, given_tag, mpi_comm_world, mpistatus, ierr)
  806. if(ierr.ne.0)then
  807. print *,'error recv_real8_3d count=',counter,'tag=',given_tag
  808. stop
  809. end if
  810. end subroutine qmpi_recv_real8_3d
  811. subroutine qmpi_recv_real4_4d(data, source, tag)
  812. implicit none
  813. real(sp) data(:,:,:,:)
  814. integer source
  815. integer, optional :: tag
  816. integer counter, given_tag
  817. given_tag=0
  818. if(present(tag)) given_tag=tag
  819. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  820. call mpi_recv(data, counter, mpi_real, source, given_tag, mpi_comm_world, mpistatus, ierr)
  821. if(ierr.ne.0)then
  822. print *,'error recv_real4_4d count=',counter,'tag=',given_tag
  823. stop
  824. end if
  825. end subroutine qmpi_recv_real4_4d
  826. subroutine qmpi_recv_real8_4d(data, source, tag)
  827. implicit none
  828. real(dp) data(:,:,:,:)
  829. integer source
  830. integer, optional :: tag
  831. integer counter, given_tag
  832. given_tag=0
  833. if(present(tag)) given_tag=tag
  834. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  835. call mpi_recv(data, counter, mpi_double_precision, source, given_tag, mpi_comm_world, mpistatus, ierr)
  836. if(ierr.ne.0)then
  837. print *,'error recv_real8_4d count=',counter,'tag=',given_tag
  838. stop
  839. end if
  840. end subroutine qmpi_recv_real8_4d
  841. subroutine qmpi_recv_logical(data, target, tag)
  842. implicit none
  843. logical data
  844. integer target, counter, given_tag
  845. integer, optional :: tag
  846. given_tag=0
  847. if(present(tag)) given_tag=tag
  848. counter=1
  849. call mpi_recv(data, counter, mpi_logical, target, given_tag, mpi_comm_world, mpistatus, ierr)
  850. if(ierr.ne.0)then
  851. print *,'error recv_logical count=',counter,'tag=',given_tag
  852. stop
  853. end if
  854. end subroutine qmpi_recv_logical
  855. subroutine qmpi_recv_string(data, target, tag)
  856. implicit none
  857. character(*) data
  858. integer target, counter, given_tag
  859. integer, optional :: tag
  860. given_tag=0
  861. if(present(tag)) given_tag=tag
  862. counter=len(data)
  863. call mpi_recv(data, counter, mpi_character, target, given_tag, mpi_comm_world, mpistatus, ierr)
  864. if(ierr.ne.0)then
  865. print *,'error recv_string count=',counter,'tag=',given_tag
  866. stop
  867. end if
  868. end subroutine qmpi_recv_string
  869. subroutine qmpi_broadcast_string(string,root)
  870. !
  871. ! send string out to all processes. if not given
  872. ! process 0 will be used as the sender - root otherwise.
  873. !
  874. implicit none
  875. character(len=*) string
  876. integer, optional :: root
  877. integer counter,boss
  878. counter=len(string)
  879. boss=0
  880. if(present(root)) then
  881. boss=root
  882. end if
  883. call mpi_bcast(string , counter, mpi_character, boss, mpi_comm_world ,ierr)
  884. end subroutine qmpi_broadcast_string
  885. subroutine qmpi_broadcast_stringarr(data,root)
  886. implicit none
  887. character(len=*) data(:)
  888. integer, optional :: root
  889. integer counter, boss
  890. counter=len(data(1))*size(data)
  891. boss=0
  892. if(present(root)) then
  893. boss=root
  894. end if
  895. call mpi_bcast(data, counter, mpi_character, boss, mpi_comm_world ,ierr)
  896. end subroutine qmpi_broadcast_stringarr
  897. subroutine qmpi_broadcast_real4(data,root)
  898. implicit none
  899. real(4) data
  900. integer, optional :: root
  901. integer counter,boss
  902. counter=1
  903. boss=0
  904. if(present(root)) boss=root
  905. call mpi_bcast(data , counter, mpi_real, boss, mpi_comm_world, ierr)
  906. end subroutine qmpi_broadcast_real4
  907. subroutine qmpi_broadcast_real8(data,root)
  908. implicit none
  909. real(8) data
  910. integer, optional :: root
  911. integer counter,boss
  912. counter=1
  913. boss=0
  914. if(present(root)) boss=root
  915. call mpi_bcast(data , counter, mpi_double_precision, boss, mpi_comm_world, ierr)
  916. end subroutine qmpi_broadcast_real8
  917. subroutine qmpi_broadcast_integer4(data,root)
  918. implicit none
  919. integer(4) data
  920. integer, optional :: root
  921. integer counter,boss
  922. counter=1
  923. boss=0
  924. if(present(root)) boss=root
  925. call mpi_bcast(data , counter, mpi_integer, boss, mpi_comm_world, ierr)
  926. end subroutine qmpi_broadcast_integer4
  927. subroutine qmpi_broadcast_integer8(data,root)
  928. implicit none
  929. integer(8) data
  930. integer, optional :: root
  931. integer counter,boss
  932. counter=1
  933. boss=0
  934. if(present(root)) boss=root
  935. call mpi_bcast(data , counter, mpi_integer8, boss, mpi_comm_world, ierr)
  936. end subroutine qmpi_broadcast_integer8
  937. subroutine qmpi_broadcast_logical(data, root)
  938. implicit none
  939. logical data
  940. integer, optional :: root
  941. integer counter,boss
  942. counter=1
  943. boss=0
  944. if(present(root)) boss=root
  945. call mpi_bcast(data , counter, mpi_logical, boss, mpi_comm_world, ierr)
  946. end subroutine qmpi_broadcast_logical
  947. subroutine qmpi_broadcast_integer4_array1d(data,root)
  948. implicit none
  949. integer(sp) data(:)
  950. integer, optional :: root
  951. integer counter,boss
  952. counter=size(data)
  953. boss=0
  954. if(present(root)) then
  955. boss=root
  956. end if
  957. call mpi_bcast(data , counter, mpi_integer, boss, mpi_comm_world ,ierr)
  958. end subroutine qmpi_broadcast_integer4_array1d
  959. subroutine qmpi_broadcast_integer8_array1d(data,root)
  960. implicit none
  961. integer(long) data(:)
  962. integer, optional :: root
  963. integer counter,boss
  964. counter=size(data)
  965. boss=0
  966. if(present(root)) then
  967. boss=root
  968. end if
  969. call mpi_bcast(data , counter, mpi_integer8, boss, mpi_comm_world ,ierr)
  970. end subroutine qmpi_broadcast_integer8_array1d
  971. subroutine qmpi_broadcast_integer4_array2d(data,root)
  972. implicit none
  973. integer(sp) data(:,:)
  974. integer, optional :: root
  975. integer counter,boss
  976. counter=size(data,1)*size(data,2)
  977. boss=0
  978. if(present(root)) then
  979. boss=root
  980. end if
  981. call mpi_bcast(data , counter, mpi_integer, boss, mpi_comm_world ,ierr)
  982. end subroutine qmpi_broadcast_integer4_array2d
  983. subroutine qmpi_broadcast_integer8_array2d(data,root)
  984. implicit none
  985. integer(long) data(:,:)
  986. integer, optional :: root
  987. integer counter,boss
  988. counter=size(data,1)*size(data,2)
  989. boss=0
  990. if(present(root)) then
  991. boss=root
  992. end if
  993. call mpi_bcast(data , counter, mpi_integer8, boss, mpi_comm_world ,ierr)
  994. end subroutine qmpi_broadcast_integer8_array2d
  995. subroutine qmpi_broadcast_real4_array1d(data,root)
  996. implicit none
  997. real(sp) data(:)
  998. integer, optional :: root
  999. integer counter, boss
  1000. counter=size(data)
  1001. boss=0
  1002. if(present(root)) then
  1003. boss=root
  1004. end if
  1005. call mpi_bcast(data , counter, mpi_real, boss, mpi_comm_world ,ierr)
  1006. end subroutine qmpi_broadcast_real4_array1d
  1007. subroutine qmpi_broadcast_real8_array1d(data,root)
  1008. implicit none
  1009. real(dp) data(:)
  1010. integer, optional :: root
  1011. integer counter, boss
  1012. counter=size(data)
  1013. boss=0
  1014. if(present(root)) then
  1015. boss=root
  1016. end if
  1017. call mpi_bcast(data , counter, mpi_double_precision, boss, mpi_comm_world ,ierr)
  1018. end subroutine qmpi_broadcast_real8_array1d
  1019. subroutine qmpi_broadcast_real4_array2d(data,root)
  1020. implicit none
  1021. real(sp) data(:,:)
  1022. integer, optional :: root
  1023. integer counter, boss
  1024. counter=size(data,1)*size(data,2)
  1025. boss=0
  1026. if(present(root)) then
  1027. boss=root
  1028. end if
  1029. call mpi_bcast(data, counter, mpi_real, boss, mpi_comm_world ,ierr)
  1030. end subroutine qmpi_broadcast_real4_array2d
  1031. subroutine qmpi_broadcast_real8_array2d(data,root)
  1032. implicit none
  1033. real(dp) data(:,:)
  1034. integer, optional :: root
  1035. integer counter, boss
  1036. counter=size(data,1)*size(data,2)
  1037. boss=0
  1038. if(present(root)) then
  1039. boss=root
  1040. end if
  1041. call mpi_bcast(data, counter, mpi_double_precision, boss, mpi_comm_world ,ierr)
  1042. end subroutine qmpi_broadcast_real8_array2d
  1043. subroutine qmpi_broadcast_real4_array3d(data,root)
  1044. implicit none
  1045. real(sp) data(:,:,:)
  1046. integer, optional :: root
  1047. integer counter, boss
  1048. counter=size(data,1)*size(data,2)*size(data,3)
  1049. boss=0
  1050. if(present(root)) then
  1051. boss=root
  1052. end if
  1053. call mpi_bcast(data , counter, mpi_real, boss, mpi_comm_world ,ierr)
  1054. end subroutine qmpi_broadcast_real4_array3d
  1055. subroutine qmpi_broadcast_real8_array3d(data,root)
  1056. implicit none
  1057. real(dp) data(:,:,:)
  1058. integer, optional :: root
  1059. integer counter, boss
  1060. counter=size(data,1)*size(data,2)*size(data,3)
  1061. boss=0
  1062. if(present(root)) then
  1063. boss=root
  1064. end if
  1065. call mpi_bcast(data , counter, mpi_double_precision, boss, mpi_comm_world ,ierr)
  1066. end subroutine qmpi_broadcast_real8_array3d
  1067. subroutine qmpi_broadcast_real4_array4d(data,root)
  1068. implicit none
  1069. real(sp) data(:,:,:,:)
  1070. integer, optional :: root
  1071. integer counter, boss
  1072. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  1073. boss=0
  1074. if(present(root)) then
  1075. boss=root
  1076. end if
  1077. call mpi_bcast(data , counter, mpi_real, boss, mpi_comm_world ,ierr)
  1078. end subroutine qmpi_broadcast_real4_array4d
  1079. subroutine qmpi_broadcast_real8_array4d(data,root)
  1080. implicit none
  1081. real(dp) data(:,:,:,:)
  1082. integer, optional :: root
  1083. integer counter, boss
  1084. counter=size(data,1)*size(data,2)*size(data,3)*size(data,4)
  1085. boss=0
  1086. if(present(root)) then
  1087. boss=root
  1088. end if
  1089. call mpi_bcast(data , counter, mpi_double_precision, boss, mpi_comm_world ,ierr)
  1090. end subroutine qmpi_broadcast_real8_array4d
  1091. subroutine qmpi_broadcast_real4s(a,b,c,d,e,f,root)
  1092. !
  1093. ! send a,b,c,d,e,f out to all processes. if not given
  1094. ! process 0 will be used as the sender - root otherwise.
  1095. !
  1096. implicit none
  1097. real(sp) a
  1098. real(sp), optional :: b,c,d,e,f
  1099. integer, optional :: root
  1100. integer counter,boss
  1101. real(sp) rbuff(6)
  1102. counter=0 ; boss=0
  1103. if(present(root)) then
  1104. boss=root
  1105. end if
  1106. ! if(present(a)) then
  1107. counter=counter+1
  1108. rbuff(counter)=a
  1109. ! end if
  1110. if(present(b)) then
  1111. counter=counter+1
  1112. rbuff(counter)=b
  1113. end if
  1114. if(present(c)) then
  1115. counter=counter+1
  1116. rbuff(counter)=c
  1117. end if
  1118. if(present(d)) then
  1119. counter=counter+1
  1120. rbuff(counter)=d
  1121. end if
  1122. if(present(e)) then
  1123. counter=counter+1
  1124. rbuff(counter)=e
  1125. end if
  1126. if(present(f)) then
  1127. counter=counter+1
  1128. rbuff(counter)=f
  1129. end if
  1130. call mpi_bcast(rbuff , counter, mpi_real, boss, mpi_comm_world ,ierr)
  1131. counter=1
  1132. a=rbuff(counter)
  1133. if(present(b)) then
  1134. counter=counter+1
  1135. b=rbuff(counter)
  1136. end if
  1137. if(present(c)) then
  1138. counter=counter+1
  1139. c=rbuff(counter)
  1140. end if
  1141. if(present(d)) then
  1142. counter=counter+1
  1143. d=rbuff(counter)
  1144. end if
  1145. if(present(e)) then
  1146. counter=counter+1
  1147. e=rbuff(counter)
  1148. end if
  1149. if(present(f)) then
  1150. counter=counter+1
  1151. f=rbuff(counter)
  1152. end if
  1153. end subroutine qmpi_broadcast_real4s
  1154. subroutine qmpi_broadcast_real8s(a,b,c,d,e,f,root)
  1155. !
  1156. ! send a,b,c,d,e,f out to all processes. if not given
  1157. ! process 0 will be used as the sender - root otherwise.
  1158. !
  1159. implicit none
  1160. real(dp) a
  1161. real(dp), optional :: b,c,d,e,f
  1162. integer, optional :: root
  1163. integer counter,boss
  1164. real(kind=8) rbuff(6)
  1165. boss=0
  1166. if(present(root)) then
  1167. boss=root
  1168. end if
  1169. counter=1
  1170. rbuff(counter)=a
  1171. if(present(b)) then
  1172. counter=counter+1
  1173. rbuff(counter)=b
  1174. end if
  1175. if(present(c)) then
  1176. counter=counter+1
  1177. rbuff(counter)=c
  1178. end if
  1179. if(present(d)) then
  1180. counter=counter+1
  1181. rbuff(counter)=d
  1182. end if
  1183. if(present(e)) then
  1184. counter=counter+1
  1185. rbuff(counter)=e
  1186. end if
  1187. if(present(f)) then
  1188. counter=counter+1
  1189. rbuff(counter)=f
  1190. end if
  1191. call mpi_bcast(rbuff , counter, mpi_double_precision, boss, mpi_comm_world ,ierr)
  1192. counter=1
  1193. a=rbuff(counter)
  1194. if(present(b)) then
  1195. counter=counter+1
  1196. b=rbuff(counter)
  1197. end if
  1198. if(present(c)) then
  1199. counter=counter+1
  1200. c=rbuff(counter)
  1201. end if
  1202. if(present(d)) then
  1203. counter=counter+1
  1204. d=rbuff(counter)
  1205. end if
  1206. if(present(e)) then
  1207. counter=counter+1
  1208. e=rbuff(counter)
  1209. end if
  1210. if(present(f)) then
  1211. counter=counter+1
  1212. f=rbuff(counter)
  1213. end if
  1214. end subroutine qmpi_broadcast_real8s
  1215. subroutine qmpi_broadcast_logicals(a,b,c,d,e,f,root)
  1216. !
  1217. ! send a,b,c,d,e,f out to all processes. if not given
  1218. ! process 0 will be used as the sender - root otherwise.
  1219. !
  1220. implicit none
  1221. logical a
  1222. logical, optional :: b,c,d,e,f
  1223. integer, optional :: root
  1224. integer counter,boss
  1225. logical lbuff(6)
  1226. boss=0
  1227. if(present(root)) then
  1228. boss=root
  1229. end if
  1230. counter=1
  1231. lbuff(counter)=a
  1232. if(present(b)) then
  1233. counter=counter+1
  1234. lbuff(counter)=b
  1235. end if
  1236. if(present(c)) then
  1237. counter=counter+1
  1238. lbuff(counter)=c
  1239. end if
  1240. if(present(d)) then
  1241. counter=counter+1
  1242. lbuff(counter)=d
  1243. end if
  1244. if(present(e)) then
  1245. counter=counter+1
  1246. lbuff(counter)=e
  1247. end if
  1248. if(present(f)) then
  1249. counter=counter+1
  1250. lbuff(counter)=f
  1251. end if
  1252. call mpi_bcast(lbuff , counter, mpi_logical, boss, mpi_comm_world ,ierr)
  1253. counter=1
  1254. a=lbuff(counter)
  1255. if(present(b)) then
  1256. counter=counter+1
  1257. b=lbuff(counter)
  1258. end if
  1259. if(present(c)) then
  1260. counter=counter+1
  1261. c=lbuff(counter)
  1262. end if
  1263. if(present(d)) then
  1264. counter=counter+1
  1265. d=lbuff(counter)
  1266. end if
  1267. if(present(e)) then
  1268. counter=counter+1
  1269. e=lbuff(counter)
  1270. end if
  1271. if(present(f)) then
  1272. counter=counter+1
  1273. f=lbuff(counter)
  1274. end if
  1275. end subroutine qmpi_broadcast_logicals
  1276. subroutine qmpi_broadcast_integer4s(a,b,c,d,e,f,root)
  1277. !
  1278. ! send a,b,c,d,e,f out to all processes. if not given
  1279. ! process 0 will be used as the sender - root otherwise.
  1280. !
  1281. implicit none
  1282. integer(sp) a
  1283. integer(sp), optional :: b,c,d,e,f,root
  1284. integer counter,boss
  1285. integer ibuff(6)
  1286. boss=0
  1287. if(present(root)) then
  1288. boss=root
  1289. end if
  1290. counter=1
  1291. ! if(present(a)) then
  1292. ! counter=counter+1
  1293. ibuff(counter)=a
  1294. ! end if
  1295. if(present(b)) then
  1296. counter=counter+1
  1297. ibuff(counter)=b
  1298. end if
  1299. if(present(c)) then
  1300. counter=counter+1
  1301. ibuff(counter)=c
  1302. end if
  1303. if(present(d)) then
  1304. counter=counter+1
  1305. ibuff(counter)=d
  1306. end if
  1307. if(present(e)) then
  1308. counter=counter+1
  1309. ibuff(counter)=e
  1310. end if
  1311. if(present(f)) then
  1312. counter=counter+1
  1313. ibuff(counter)=f
  1314. end if
  1315. call mpi_bcast(ibuff , counter, mpi_integer, boss, mpi_comm_world ,ierr)
  1316. counter=1
  1317. a=ibuff(counter)
  1318. if(present(b)) then
  1319. counter=counter+1
  1320. b=ibuff(counter)
  1321. end if
  1322. if(present(c)) then
  1323. counter=counter+1
  1324. c=ibuff(counter)
  1325. end if
  1326. if(present(d)) then
  1327. counter=counter+1
  1328. d=ibuff(counter)
  1329. end if
  1330. if(present(e)) then
  1331. counter=counter+1
  1332. e=ibuff(counter)
  1333. end if
  1334. if(present(f)) then
  1335. counter=counter+1
  1336. f=ibuff(counter)
  1337. end if
  1338. end subroutine qmpi_broadcast_integer4s
  1339. subroutine qmpi_broadcast_integer8s(a,b,c,d,e,f,root)
  1340. !
  1341. ! send a,b,c,d,e,f out to all processes. if not given
  1342. ! process 0 will be used as the sender - root otherwise.
  1343. !
  1344. implicit none
  1345. integer(long) a
  1346. integer(long), optional :: b,c,d,e,f,root
  1347. integer counter,boss
  1348. integer ibuff(6)
  1349. boss=0
  1350. if(present(root)) then
  1351. boss=root
  1352. end if
  1353. counter=1
  1354. ! if(present(a)) then
  1355. ! counter=counter+1
  1356. ibuff(counter)=a
  1357. ! end if
  1358. if(present(b)) then
  1359. counter=counter+1
  1360. ibuff(counter)=b
  1361. end if
  1362. if(present(c)) then
  1363. counter=counter+1
  1364. ibuff(counter)=c
  1365. end if
  1366. if(present(d)) then
  1367. counter=counter+1
  1368. ibuff(counter)=d
  1369. end if
  1370. if(present(e)) then
  1371. counter=counter+1
  1372. ibuff(counter)=e
  1373. end if
  1374. if(present(f)) then
  1375. counter=counter+1
  1376. ibuff(counter)=f
  1377. end if
  1378. call mpi_bcast(ibuff , counter, mpi_integer8, boss, mpi_comm_world ,ierr)
  1379. counter=1
  1380. a=ibuff(counter)
  1381. if(present(b)) then
  1382. counter=counter+1
  1383. b=ibuff(counter)
  1384. end if
  1385. if(present(c)) then
  1386. counter=counter+1
  1387. c=ibuff(counter)
  1388. end if
  1389. if(present(d)) then
  1390. counter=counter+1
  1391. d=ibuff(counter)
  1392. end if
  1393. if(present(e)) then
  1394. counter=counter+1
  1395. e=ibuff(counter)
  1396. end if
  1397. if(present(f)) then
  1398. counter=counter+1
  1399. f=ibuff(counter)
  1400. end if
  1401. end subroutine qmpi_broadcast_integer8s
  1402. subroutine qmpi_real_reduction(type,a,b,c,d,e,f,root)
  1403. !
  1404. ! perform a reduction of 'type' on each of the given arguments a - f.
  1405. ! if type is:
  1406. ! 'sum': for each argument, return the sum of the argument over all processors
  1407. ! 'mul': the product
  1408. ! 'min': the minimum value
  1409. ! 'max': the maximum value
  1410. ! root is an optional argument, if given only return the result on that processor (reduce)
  1411. ! the default is to return the result on all processors (allreduce)
  1412. !
  1413. implicit none
  1414. character(3) type
  1415. real(sp) a
  1416. real(sp), optional, intent(inout) :: b,c,d,e,f
  1417. integer, optional :: root
  1418. integer counter,boss
  1419. integer, parameter :: dp=8
  1420. real(dp) rbuff(6),globrbuff(6)
  1421. if( trim(type).ne.'sum' .and. trim(type).ne.'mul' .and. trim(type).ne.'min' .and. trim(type).ne.'max')then
  1422. print *,'qmpi.f90 reduce error: reduction of type ',type,'not supported'
  1423. stop
  1424. end if
  1425. boss=0
  1426. if(present(root)) boss=root
  1427. globrbuff(:)=0.0
  1428. counter=0
  1429. ! if(present(a)) then
  1430. counter=counter+1
  1431. rbuff(counter)=real(a,dp)
  1432. ! end if
  1433. if(present(b)) then
  1434. counter=counter+1
  1435. rbuff(counter)=real(b,dp)
  1436. end if
  1437. if(present(c)) then
  1438. counter=counter+1
  1439. rbuff(counter)=real(c,dp)
  1440. end if
  1441. if(present(d)) then
  1442. counter=counter+1
  1443. rbuff(counter)=real(d,dp)
  1444. end if
  1445. if(present(e)) then
  1446. counter=counter+1
  1447. rbuff(counter)=real(e,dp)
  1448. end if
  1449. if(present(f)) then
  1450. counter=counter+1
  1451. rbuff(counter)=real(f,dp)
  1452. end if
  1453. select case(type)
  1454. case('sum')
  1455. if(present(root))then
  1456. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_sum,boss,mpi_comm_world,ierr)
  1457. else
  1458. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_sum,mpi_comm_world,ierr)
  1459. end if
  1460. case('mul')
  1461. if(present(root))then
  1462. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_prod,boss,mpi_comm_world,ierr)
  1463. else
  1464. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_prod,mpi_comm_world,ierr)
  1465. end if
  1466. case('min')
  1467. if(present(root))then
  1468. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_min,boss,mpi_comm_world,ierr)
  1469. else
  1470. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_min,mpi_comm_world,ierr)
  1471. end if
  1472. case('max')
  1473. if(present(root))then
  1474. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_max,boss,mpi_comm_world,ierr)
  1475. else
  1476. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_max,mpi_comm_world,ierr)
  1477. end if
  1478. end select
  1479. counter=0
  1480. ! if(present(a)) then
  1481. counter=counter+1
  1482. a=globrbuff(counter)
  1483. ! end if
  1484. if(present(b)) then
  1485. counter=counter+1
  1486. b=globrbuff(counter)
  1487. end if
  1488. if(present(c)) then
  1489. counter=counter+1
  1490. c=globrbuff(counter)
  1491. end if
  1492. if(present(d)) then
  1493. counter=counter+1
  1494. d=globrbuff(counter)
  1495. end if
  1496. if(present(e)) then
  1497. counter=counter+1
  1498. e=globrbuff(counter)
  1499. end if
  1500. if(present(f)) then
  1501. counter=counter+1
  1502. f=globrbuff(counter)
  1503. end if
  1504. end subroutine qmpi_real_reduction
  1505. subroutine qmpi_real8_reduction(type,a,b,c,d,e,f,root)
  1506. !
  1507. ! perform a reduction of 'type' on each of the given arguments a - f.
  1508. ! if type is:
  1509. ! 'sum': for each argument, return the sum of the argument over all processors
  1510. ! 'mul': the product
  1511. ! 'min': the minimum value
  1512. ! 'max': the maximum value
  1513. ! root is an optional argument, if given only return the result on that processor (reduce)
  1514. ! the default is to return the result on all processors (allreduce)
  1515. !
  1516. implicit none
  1517. integer, parameter :: dp=8
  1518. character(3) type
  1519. real(dp) a
  1520. real(dp), optional, intent(inout) :: b,c,d,e,f
  1521. integer, optional :: root
  1522. integer counter,boss
  1523. real(dp) rbuff(6),globrbuff(6)
  1524. if( trim(type).ne.'sum' .and. trim(type).ne.'mul' .and. trim(type).ne.'min' .and. trim(type).ne.'max')then
  1525. print *,'qmpi.f90 reduce error: reduction of type ',type,'not supported'
  1526. stop
  1527. end if
  1528. boss=0
  1529. if(present(root))boss=root
  1530. globrbuff(:)=0.0
  1531. counter=0
  1532. ! if(present(a)) then
  1533. counter=counter+1
  1534. rbuff(counter)=a
  1535. ! end if
  1536. if(present(b)) then
  1537. counter=counter+1
  1538. rbuff(counter)=b
  1539. end if
  1540. if(present(c)) then
  1541. counter=counter+1
  1542. rbuff(counter)=c
  1543. end if
  1544. if(present(d)) then
  1545. counter=counter+1
  1546. rbuff(counter)=d
  1547. end if
  1548. if(present(e)) then
  1549. counter=counter+1
  1550. rbuff(counter)=e
  1551. end if
  1552. if(present(f)) then
  1553. counter=counter+1
  1554. rbuff(counter)=f
  1555. end if
  1556. select case(type)
  1557. case('sum')
  1558. if(present(root))then
  1559. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_sum,boss,mpi_comm_world,ierr)
  1560. else
  1561. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_sum,mpi_comm_world,ierr)
  1562. end if
  1563. case('mul')
  1564. if(present(root))then
  1565. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_prod,boss,mpi_comm_world,ierr)
  1566. else
  1567. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_prod,mpi_comm_world,ierr)
  1568. end if
  1569. case('min')
  1570. if(present(root))then
  1571. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_min,boss,mpi_comm_world,ierr)
  1572. else
  1573. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_min,mpi_comm_world,ierr)
  1574. end if
  1575. case('max')
  1576. if(present(root))then
  1577. call mpi_reduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_max,boss,mpi_comm_world,ierr)
  1578. else
  1579. call mpi_allreduce(rbuff,globrbuff,counter,mpi_double_precision,mpi_max,mpi_comm_world,ierr)
  1580. end if
  1581. end select
  1582. counter=0
  1583. ! if(present(a)) then
  1584. counter=counter+1
  1585. a=globrbuff(counter)
  1586. ! end if
  1587. if(present(b)) then
  1588. counter=counter+1
  1589. b=globrbuff(counter)
  1590. end if
  1591. if(present(c)) then
  1592. counter=counter+1
  1593. c=globrbuff(counter)
  1594. end if
  1595. if(present(d)) then
  1596. counter=counter+1
  1597. d=globrbuff(counter)
  1598. end if
  1599. if(present(e)) then
  1600. counter=counter+1
  1601. e=globrbuff(counter)
  1602. end if
  1603. if(present(f)) then
  1604. counter=counter+1
  1605. f=globrbuff(counter)
  1606. end if
  1607. end subroutine qmpi_real8_reduction
  1608. subroutine qmpi_integer_reduction(type,a,b,c,d,e,f,root)
  1609. !
  1610. ! perform a reduction of 'type' on each of the given arguments a - f.
  1611. ! if type is:
  1612. ! 'sum': for each argument, return the sum of the argument over all processors
  1613. ! 'mul': the product
  1614. ! 'min': the minimum value
  1615. ! 'max': the maximum value
  1616. ! root is an optional argument, if given only return the result on that processor (reduce)
  1617. ! the default is to return the result on all processors (allreduce)
  1618. !
  1619. implicit none
  1620. character(3) type
  1621. integer(sp) a
  1622. integer(sp), optional, intent(inout) :: b,c,d,e,f
  1623. integer, optional :: root
  1624. integer counter,boss
  1625. integer rbuff(6),globrbuff(6)
  1626. if( trim(type).ne.'sum' .and. trim(type).ne.'mul' .and. trim(type).ne.'min' .and. trim(type).ne.'max')then
  1627. print *,'qmpi.f90 reduce error: reduction of type ',type,'not supported'
  1628. stop
  1629. end if
  1630. boss=0
  1631. if(present(root))boss=root
  1632. globrbuff(:)=0
  1633. counter=0
  1634. !if(present(a)) then
  1635. counter=counter+1
  1636. rbuff(counter)=a
  1637. !end if
  1638. if(present(b)) then
  1639. counter=counter+1
  1640. rbuff(counter)=b
  1641. end if
  1642. if(present(c)) then
  1643. counter=counter+1
  1644. rbuff(counter)=c
  1645. end if
  1646. if(present(d)) then
  1647. counter=counter+1
  1648. rbuff(counter)=d
  1649. end if
  1650. if(present(e)) then
  1651. counter=counter+1
  1652. rbuff(counter)=e
  1653. end if
  1654. if(present(f)) then
  1655. counter=counter+1
  1656. rbuff(counter)=f
  1657. end if
  1658. select case(type)
  1659. case('sum')
  1660. if(present(root))then
  1661. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_sum,boss,mpi_comm_world,ierr)
  1662. else
  1663. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_sum,mpi_comm_world,ierr)
  1664. end if
  1665. case('mul')
  1666. if(present(root))then
  1667. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_prod,boss,mpi_comm_world,ierr)
  1668. else
  1669. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_prod,mpi_comm_world,ierr)
  1670. end if
  1671. case('min')
  1672. if(present(root))then
  1673. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_min,boss,mpi_comm_world,ierr)
  1674. else
  1675. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_min,mpi_comm_world,ierr)
  1676. end if
  1677. case('max')
  1678. if(present(root))then
  1679. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_max,boss,mpi_comm_world,ierr)
  1680. else
  1681. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER, mpi_max,mpi_comm_world,ierr)
  1682. end if
  1683. end select
  1684. counter=0
  1685. ! if(present(a)) then
  1686. counter=counter+1
  1687. a=globrbuff(counter)
  1688. ! end if
  1689. if(present(b)) then
  1690. counter=counter+1
  1691. b=globrbuff(counter)
  1692. end if
  1693. if(present(c)) then
  1694. counter=counter+1
  1695. c=globrbuff(counter)
  1696. end if
  1697. if(present(d)) then
  1698. counter=counter+1
  1699. d=globrbuff(counter)
  1700. end if
  1701. if(present(e)) then
  1702. counter=counter+1
  1703. e=globrbuff(counter)
  1704. end if
  1705. if(present(f)) then
  1706. counter=counter+1
  1707. f=globrbuff(counter)
  1708. end if
  1709. end subroutine qmpi_integer_reduction
  1710. subroutine qmpi_integer8_reduction(type,a,b,c,d,e,f,root)
  1711. !
  1712. ! perform a reduction of 'type' on each of the given arguments a - f.
  1713. ! if type is:
  1714. ! 'sum': for each argument, return the sum of the argument over all processors
  1715. ! 'mul': the product
  1716. ! 'min': the minimum value
  1717. ! 'max': the maximum value
  1718. ! root is an optional argument, if given only return the result on that processor (reduce)
  1719. ! the default is to return the result on all processors (allreduce)
  1720. !
  1721. implicit none
  1722. character(3) type
  1723. integer(long) a
  1724. integer(long), optional, intent(inout) :: b,c,d,e,f
  1725. integer, optional :: root
  1726. integer counter,boss
  1727. integer(long) rbuff(6),globrbuff(6)
  1728. if(len(type).ne.3)then
  1729. print *,'qmpi.f90 reduce error: type must be one of "mul","sum","min" or "max"'
  1730. stop
  1731. end if
  1732. if( trim(type).ne.'sum' .and. trim(type).ne.'mul' .and. trim(type).ne.'min' .and. trim(type).ne.'max')then
  1733. print *,'qmpi.f90 reduce error: reduction of type ',type,'not supported'
  1734. stop
  1735. end if
  1736. boss=0
  1737. if(present(root))boss=root
  1738. globrbuff(:)=0_dp
  1739. counter=0
  1740. ! if(present(a)) then
  1741. counter=counter+1
  1742. rbuff(counter)=a
  1743. ! end if
  1744. if(present(b)) then
  1745. counter=counter+1
  1746. rbuff(counter)=b
  1747. end if
  1748. if(present(c)) then
  1749. counter=counter+1
  1750. rbuff(counter)=c
  1751. end if
  1752. if(present(d)) then
  1753. counter=counter+1
  1754. rbuff(counter)=d
  1755. end if
  1756. if(present(e)) then
  1757. counter=counter+1
  1758. rbuff(counter)=e
  1759. end if
  1760. if(present(f)) then
  1761. counter=counter+1
  1762. rbuff(counter)=f
  1763. end if
  1764. select case(type)
  1765. case('sum')
  1766. if(present(root))then
  1767. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_sum,boss,mpi_comm_world,ierr)
  1768. else
  1769. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_sum,mpi_comm_world,ierr)
  1770. end if
  1771. case('mul')
  1772. if(present(root))then
  1773. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_prod,boss,mpi_comm_world,ierr)
  1774. else
  1775. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_prod,mpi_comm_world,ierr)
  1776. end if
  1777. case('min')
  1778. if(present(root))then
  1779. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_min,boss,mpi_comm_world,ierr)
  1780. else
  1781. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_min,mpi_comm_world,ierr)
  1782. end if
  1783. case('max')
  1784. if(present(root))then
  1785. call mpi_reduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_max,boss,mpi_comm_world,ierr)
  1786. else
  1787. call mpi_allreduce(rbuff,globrbuff,counter, MPI_INTEGER8, mpi_max,mpi_comm_world,ierr)
  1788. end if
  1789. end select
  1790. counter=1
  1791. a=globrbuff(counter)
  1792. if(present(b)) then
  1793. counter=counter+1
  1794. b=globrbuff(counter)
  1795. end if
  1796. if(present(c)) then
  1797. counter=counter+1
  1798. c=globrbuff(counter)
  1799. end if
  1800. if(present(d)) then
  1801. counter=counter+1
  1802. d=globrbuff(counter)
  1803. end if
  1804. if(present(e)) then
  1805. counter=counter+1
  1806. e=globrbuff(counter)
  1807. end if
  1808. if(present(f)) then
  1809. counter=counter+1
  1810. f=globrbuff(counter)
  1811. end if
  1812. end subroutine qmpi_integer8_reduction
  1813. ! later?
  1814. ! packing to reduce number of sends:
  1815. ! call pack(u)
  1816. ! call pack(eta(1,:))
  1817. ! call pack(v)
  1818. ! call send_pack(1)
  1819. ! ...
  1820. ! call receive_pack(0)
  1821. ! call unpack(u)
  1822. ! call unpack(eta(1,:)
  1823. !
  1824. end module qmpi
  1825. #else
  1826. #warning "COMPILING WITHOUT QMPI CODE"
  1827. module qmpi_fake
  1828. implicit none
  1829. logical, parameter :: master = .true.
  1830. integer, parameter :: qmpi_num_proc = 1
  1831. integer, parameter :: qmpi_proc_num = 0
  1832. contains
  1833. subroutine stop_mpi()
  1834. stop
  1835. end subroutine stop_mpi
  1836. end module qmpi_fake
  1837. #endif