Upped the version to 3.2.0
[gromacs.git] / src / gmxlib / mvxvf.c
blob042bc607becaff851ab8a9b7534795af5ad24e6a
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 /* This file is completely threadsafe - keep it that way! */
38 #include <sysstuff.h>
39 #include <string.h>
40 #include "typedefs.h"
41 #include "main.h"
42 #include "assert.h"
43 #include "mvdata.h"
44 #include "network.h"
45 #include "smalloc.h"
46 #include "fatal.h"
47 #include "symtab.h"
48 #include "main.h"
49 #include "typedefs.h"
50 #include "vec.h"
51 #include "tgroup.h"
52 #include "block_tx.h"
53 #include "nrnb.h"
55 void move_rvecs(FILE *log,bool bForward,bool bSum,
56 int left,int right,rvec vecs[],rvec buf[],
57 int shift,t_nsborder *nsb,t_nrnb *nrnb)
59 int i,j,j0=137,j1=391;
60 int cur,nsum;
61 #define next ((cur+1) % nsb->nnodes)
62 #define prev ((cur-1+nsb->nnodes) % nsb->nnodes)
64 if (bSum)
65 cur=(nsb->nodeid+nsb->shift) % nsb->nnodes;
66 else
67 cur=nsb->nodeid;
69 nsum=0;
70 for(i=0; (i<shift); i++) {
71 if (bSum) {
72 if (bForward) {
73 j0=nsb->index[prev];
74 j1=j0+nsb->homenr[prev];
76 else {
77 j0=nsb->index[next];
78 j1=j0+nsb->homenr[next];
80 for(j=j0; (j<j1); j++) {
81 clear_rvec(buf[j]);
84 /* Forward pulse around the ring, to increasing NODE number */
85 if (bForward) {
86 if (bSum)
87 gmx_tx_rx_real(right,vecs[nsb->index[cur]], nsb->homenr[cur]*DIM,
88 left,buf [nsb->index[prev]],nsb->homenr[prev]*DIM);
89 else
90 gmx_tx_rx_real(right,vecs[nsb->index[cur]], nsb->homenr[cur]*DIM,
91 left, vecs[nsb->index[prev]],nsb->homenr[prev]*DIM);
92 /* Wait for communication to end */
93 gmx_wait(right,left);
96 /* Backward pulse around the ring, to decreasing NODE number */
97 else {
98 if (bSum)
99 gmx_tx_rx_real(left, vecs[nsb->index[cur]], nsb->homenr[cur]*DIM,
100 right,buf [nsb->index[next]],nsb->homenr[next]*DIM);
101 else
102 gmx_tx_rx_real(left, vecs[nsb->index[cur]], nsb->homenr[cur]*DIM,
103 right,vecs[nsb->index[next]],nsb->homenr[next]*DIM);
104 /* Wait for communication to end */
105 gmx_wait(left,right);
108 /* Actual summation */
109 if (bSum) {
110 for(j=j0; (j<j1); j++) {
111 rvec_inc(vecs[j],buf[j]);
113 nsum+=(j1-j0);
115 if (bForward)
116 cur=prev;
117 else
118 cur=next;
120 if (nsum > 0)
121 inc_nrnb(nrnb,eNR_FSUM,nsum);
122 #undef next
123 #undef prev
126 void move_x(FILE *log,
127 int left,int right,rvec x[],t_nsborder *nsb,
128 t_nrnb *nrnb)
130 move_rvecs(log,FALSE,FALSE,left,right,x,NULL,nsb->shift,nsb,nrnb);
131 move_rvecs(log,TRUE, FALSE,left,right,x,NULL,nsb->bshift,nsb,nrnb);
133 where();
136 void move_f(FILE *log,
137 int left,int right,rvec f[],rvec fadd[],
138 t_nsborder *nsb,t_nrnb *nrnb)
140 move_rvecs(log,TRUE, TRUE,left,right,f,fadd,nsb->shift,nsb,nrnb);
141 move_rvecs(log,FALSE,TRUE,left,right,f,fadd,nsb->bshift,nsb,nrnb);
143 where();
146 void move_cgcm(FILE *log,t_commrec *cr,rvec cg_cm[],int nload[])
148 int i,start,nr;
149 int cur=cr->nodeid;
150 #define next ((cur+1) % cr->nnodes)
152 for(i=0; (i<cr->nnodes-1); i++) {
153 start = (cur == 0) ? 0 : nload[cur-1];
154 nr = nload[cur] - start;
155 gmx_tx(cr->left, cg_cm[start], nr*sizeof(cg_cm[0]));
156 #ifdef DEBUG
157 fprintf(log,"move_cgcm: TX start=%d, nr=%d\n",start,nr);
158 #endif
159 start = (next == 0) ? 0 : nload[next-1];
160 nr = nload[next] - start;
161 gmx_rx(cr->right,cg_cm[start], nr*sizeof(cg_cm[0]));
162 #ifdef DEBUG
163 fprintf(log,"move_cgcm: RX start=%d, nr=%d\n",start,nr);
164 #endif
165 gmx_tx_wait(cr->left);
166 gmx_rx_wait(cr->right);
168 cur=next;
170 #undef next