exchange_s.c 5.95 KB
Newer Older
Tilman Steinweg's avatar
Tilman Steinweg committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41
/*-----------------------------------------------------------------------------------------
 * Copyright (C) 2013  For the list of authors, see file AUTHORS.
 *
 * This file is part of DENISE.
 * 
 * DENISE is free software: you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 * 
 * DENISE is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with DENISE. See file COPYING and/or <http://www.gnu.org/licenses/gpl-2.0.html>.
-----------------------------------------------------------------------------------------*/

/*------------------------------------------------------------------------
 *   write values of dynamic field variables at the edges of the
 *   local grid into buffer arrays and  exchanged between
 *   processes.
 *   last update 21/09/02, T. Bohlen
 *
 *  ----------------------------------------------------------------------*/

#include "fd.h"

void exchange_s(float ** sxx, float ** syy, 
	float ** sxy, float ** bufferlef_to_rig, float ** bufferrig_to_lef, 
	float ** buffertop_to_bot, float ** bufferbot_to_top,
	MPI_Request * req_send, MPI_Request * req_rec){


	extern int NX, NY, POS[3], NPROCX, NPROCY, BOUNDARY, FDORDER;
	extern int INDEX[5];
	extern const int TAG1,TAG2,TAG5,TAG6;
	MPI_Status  status;
	int i, j, fdo, fdo3, n, l;


42
	fdo = FDORDER/2;
Tilman Steinweg's avatar
Tilman Steinweg committed
43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127
	fdo3 = 2*fdo;
	
	
	/* top - bottom */
	
	if (POS[2]!=0)	/* no boundary exchange at top of global grid */
	for (i=1;i<=NX;i++){
		/* storage of top of local volume into buffer */
		n = 1;
		for (l=1;l<=fdo-1;l++) {
			buffertop_to_bot[i][n++]  = sxy[l][i];
		}
		for (l=1;l<=fdo;l++) {
			buffertop_to_bot[i][n++]  = syy[l][i];
		}
	}


	if (POS[2]!=NPROCY-1)	/* no boundary exchange at bottom of global grid */
	for (i=1;i<=NX;i++){			
		/* storage of bottom of local volume into buffer */
		n = 1;
		for (l=1;l<=fdo;l++) {
			bufferbot_to_top[i][n++]  = sxy[NY-l+1][i];
		}
		for (l=1;l<=fdo-1;l++) {
			bufferbot_to_top[i][n++]  = syy[NY-l+1][i];
		}
	}
	
	
  	 /* send and receive values for points at inner boundaries */

/*
	MPI_Bsend(&buffertop_to_bot[1][1],NX*fdo3,MPI_FLOAT,INDEX[3],TAG5,MPI_COMM_WORLD);
	MPI_Barrier(MPI_COMM_WORLD);
	MPI_Recv(&buffertop_to_bot[1][1],NX*fdo3,MPI_FLOAT,INDEX[4],TAG5,MPI_COMM_WORLD,&status);
	MPI_Bsend(&bufferbot_to_top[1][1],NX*fdo3,MPI_FLOAT,INDEX[4],TAG6,MPI_COMM_WORLD);
	MPI_Barrier(MPI_COMM_WORLD);
	MPI_Recv(&bufferbot_to_top[1][1],NX*fdo3,MPI_FLOAT,INDEX[3],TAG6,MPI_COMM_WORLD,&status);   
*/
	/* Initiates a communication with a persistent request handle */
	/*for (i=2;i<=3;i++){
		MPI_Start(&req_send[i]);
		MPI_Wait(&req_send[i],&status);
		MPI_Start(&req_rec[i]);
		MPI_Wait(&req_rec[i],&status);
	}*/
	
	/* alternative communication */
	/* still blocking communication */
	MPI_Sendrecv_replace(&buffertop_to_bot[1][1],NX*fdo3,MPI_FLOAT,INDEX[3],TAG5,INDEX[4],TAG5,MPI_COMM_WORLD,&status);
	MPI_Sendrecv_replace(&bufferbot_to_top[1][1],NX*fdo3,MPI_FLOAT,INDEX[4],TAG6,INDEX[3],TAG6,MPI_COMM_WORLD,&status);
	


	if (POS[2]!=NPROCY-1)	/* no boundary exchange at bottom of global grid */
	for (i=1;i<=NX;i++){
		n = 1;
		for (l=1;l<=fdo-1;l++) {
			sxy[NY+l][i] = buffertop_to_bot[i][n++];
		}
		for (l=1;l<=fdo;l++) {
			syy[NY+l][i] = buffertop_to_bot[i][n++];
		}
	}

	if (POS[2]!=0)	/* no boundary exchange at top of global grid */
	for (i=1;i<=NX;i++){
		n = 1;
		for (l=1;l<=fdo;l++) {
			sxy[1-l][i] = bufferbot_to_top[i][n++];
		}
		for (l=1;l<=fdo-1;l++) {
			syy[1-l][i] = bufferbot_to_top[i][n++];
		}
	}
	
	
	/* left - right */

	if ((BOUNDARY) || (POS[1]!=0))	/* no boundary exchange at left edge of global grid */
	for (j=1;j<=NY;j++){
		/* storage of left edge of local volume into buffer */
		n = 1;
128
		for (l=1;l<=fdo-1;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
129 130
			bufferlef_to_rig[j][n++] =  sxy[j][l];
		}
131
		for (l=1;l<=fdo;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
132 133 134 135 136 137 138 139 140
			bufferlef_to_rig[j][n++] =  sxx[j][l];
		}
	}


	if ((BOUNDARY) || (POS[1]!=NPROCX-1))	/* no boundary exchange at right edge of global grid */
	for (j=1;j<=NY;j++){
		/* storage of right edge of local volume into buffer */
		n = 1;
141
		for (l=1;l<=fdo;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
142 143
			bufferrig_to_lef[j][n++] =  sxy[j][NX-l+1];
		}
144
		for (l=1;l<=fdo-1;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180
			bufferrig_to_lef[j][n++] =  sxx[j][NX-l+1];
		}
	}	
	
	


 	 /* send and receive values for points at inner boundaries */

/*
 	MPI_Bsend(&bufferlef_to_rig[1][1],(NY)*fdo3,MPI_FLOAT,INDEX[1],TAG1,MPI_COMM_WORLD);
	MPI_Barrier(MPI_COMM_WORLD);
	MPI_Recv(&bufferlef_to_rig[1][1],(NY)*fdo3,MPI_FLOAT,INDEX[2],TAG1,MPI_COMM_WORLD,&status);
	MPI_Bsend(&bufferrig_to_lef[1][1],(NY)*fdo3,MPI_FLOAT,INDEX[2],TAG2,MPI_COMM_WORLD);
	MPI_Barrier(MPI_COMM_WORLD);
	MPI_Recv(&bufferrig_to_lef[1][1],(NY)*fdo3,MPI_FLOAT,INDEX[1],TAG2,MPI_COMM_WORLD,&status);
*/

	/* send and reveive values at edges of the local grid */
	/*for (i=0;i<=1;i++){
		MPI_Start(&req_send[i]);
		MPI_Wait(&req_send[i],&status);
		MPI_Start(&req_rec[i]);
		MPI_Wait(&req_rec[i],&status);
	}*/
	
	
	/* alternative communication */
	/* still blocking communication */
	MPI_Sendrecv_replace(&bufferlef_to_rig[1][1],NY*fdo3,MPI_FLOAT,INDEX[1],TAG1,INDEX[2],TAG1,MPI_COMM_WORLD,&status);
	MPI_Sendrecv_replace(&bufferrig_to_lef[1][1],NY*fdo3,MPI_FLOAT,INDEX[2],TAG2,INDEX[1],TAG2,MPI_COMM_WORLD,&status);	

	
	if ((BOUNDARY) || (POS[1]!=NPROCX-1))	/* no boundary exchange at right edge of global grid */
	for (j=1;j<=NY;j++){
		n = 1;
181
		for (l=1;l<=fdo-1;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
182 183
			sxy[j][NX+l] = bufferlef_to_rig[j][n++];
		}
184
		for (l=1;l<=fdo;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
185 186 187 188 189 190 191
			sxx[j][NX+l] = bufferlef_to_rig[j][n++];
		}
	}

	if ((BOUNDARY) || (POS[1]!=0))	/* no boundary exchange at left edge of global grid */
	for (j=1;j<=NY;j++){
		n = 1;
192
		for (l=1;l<=fdo;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
193 194
			sxy[j][1-l] = bufferrig_to_lef[j][n++];
		}
195
		for (l=1;l<=fdo-1;l++) {
Tilman Steinweg's avatar
Tilman Steinweg committed
196 197 198 199 200 201
			sxx[j][1-l] = bufferrig_to_lef[j][n++];
		}
	}


}