source: CIVL/examples/omp/HydroC/SplitSurface.c@ 1aaefd4

main test-branch
Last change on this file since 1aaefd4 was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 7.6 KB
Line 
1/*
2 A simple 2D hydro code
3 (C) Romain Teyssier : CEA/IRFU -- original F90 code
4 (C) Pierre-Francois Lavallee : IDRIS -- original F90 code
5 (C) Guillaume Colin de Verdiere : CEA/DAM -- for the C version
6 (C) Adèle Villiermet : CINES -- for FTI integration
7*/
8
9/*
10
11This software is governed by the CeCILL license under French law and
12abiding by the rules of distribution of free software. You can use,
13modify and/ or redistribute the software under the terms of the CeCILL
14license as circulated by CEA, CNRS and INRIA at the following URL
15"http://www.cecill.info".
16
17As a counterpart to the access to the source code and rights to copy,
18modify and redistribute granted by the license, users are provided only
19with a limited warranty and the software's author, the holder of the
20economic rights, and the successive licensors have only limited
21liability.
22
23In this respect, the user's attention is drawn to the risks associated
24with loading, using, modifying and/or developing or reproducing the
25software by the user in light of its specific status of free software,
26that may mean that it is complicated to manipulate, and that also
27therefore means that it is reserved for developers and experienced
28professionals having in-depth computer knowledge. Users are therefore
29encouraged to load and test the software's suitability as regards their
30requirements in conditions enabling the security of their systems and/or
31data to be ensured and, more generally, to use and operate it in the
32same conditions as regards security.
33
34The fact that you are presently reading this means that you have had
35knowledge of the CeCILL license and that you accept its terms.
36
37*/
38#ifdef MPI
39#if FTI>0
40#include <fti.h>
41#endif
42#include <mpi.h>
43#endif
44#include <stdio.h>
45#include <string.h>
46#include <strings.h>
47#include <unistd.h>
48#include <stdlib.h>
49#include <malloc.h>
50#include <assert.h>
51#include <limits.h>
52#include <sys/time.h>
53#include <float.h>
54#include <math.h>
55
56//
57#include "SplitSurface.h"
58
59/*
60 This module splits a surface according to a k-d tree
61*/
62
63static char
64SplitSurface(int level,
65 int xmin, int xmax,
66 int ymin, int ymax,
67 int *xmin1, int *xmax1, int *ymin1, int *ymax1, int *xmin2, int *xmax2, int *ymin2, int *ymax2) {
68 char split = 0;
69
70 *xmin1 = *xmin2 = xmin;
71 *ymin1 = *ymin2 = ymin;
72 *xmax1 = *xmax2 = xmax;
73 *ymax1 = *ymax2 = ymax;
74
75 switch (level % 2) {
76 case 0:
77 if (xmin != (xmax - 1)) {
78 *xmin2 = (xmin + xmax + (level % 2)) / 2;
79 *xmax1 = *xmin2;
80 split = 'X';
81 break;
82 }
83 case 1:
84 if (ymin != (ymax - 1)) {
85 *ymin2 = (ymin + ymax + (level % 2)) / 2;
86 *ymax1 = *ymin2;
87 split = 'Y';
88 break;
89 }
90 } // switch
91
92 if (!split && ((xmin != (xmax - 1)) || (ymin != (ymax - 1)))) {
93 split = SplitSurface(level - 1, xmin, xmax, ymin, ymax, xmin1, xmax1, ymin1, ymax1, xmin2, xmax2, ymin2, ymax2);
94 }
95 return split;
96}
97
98#ifdef WITHKDTREE
99void
100CalcSubSurface(int xmin, int xmax,
101 int ymin, int ymax, int pmin, int pmax, int level, int box[MAXBOX_D], int mype, int pass) {
102 long pmid;
103 int split;
104 int xmin1, xmax1, ymin1, ymax1;
105 int xmin2, xmax2, ymin2, ymax2;
106
107 // if (mype == 0) printf("Evaluating a KDTree\n");
108
109 if (pmin == pmax) {
110 // We're down to one PE it's our sub vol limits.
111 // printf("[%4d] %4d %4d %4d %4d (%d)\n", pmin, xmin, xmax, ymin, ymax, level);
112 if (pmin == mype) {
113 box[0] = xmin;
114 box[1] = xmax;
115 box[2] = ymin;
116 box[3] = ymax;
117 } else {
118 // look for a common edge
119 if ((box[XMIN_D] == xmax) && (box[YMIN_D] == ymin) && (box[YMAX_D] == ymax))
120 box[LEFT_D] = pmin;
121 if ((box[XMAX_D] == xmin) && (box[YMIN_D] == ymin) && (box[YMAX_D] == ymax))
122 box[RIGHT_D] = pmin;
123 if ((box[YMIN_D] == ymax) && (box[XMIN_D] == xmin) && (box[XMAX_D] == xmax))
124 box[DOWN_D] = pmin;
125 if ((box[YMAX_D] == ymin) && (box[XMIN_D] == xmin) && (box[XMAX_D] == xmax))
126 box[UP_D] = pmin;
127 }
128 return;
129 }
130
131 split = SplitSurface(level, xmin, xmax, ymin, ymax, &xmin1, &xmax1, &ymin1, &ymax1, &xmin2, &xmax2, &ymin2, &ymax2);
132
133 if (split) {
134 // recurse on sub problems
135 pmid = (pmax + pmin) / 2;
136 CalcSubSurface(xmin1, xmax1, ymin1, ymax1, pmin, pmid, level + 1, box, mype, pass);
137 CalcSubSurface(xmin2, xmax2, ymin2, ymax2, pmid + 1, pmax, level + 1, box, mype, pass);
138 } else {
139 // cerr << "Too many PE for this problem size " << endl;
140 fprintf(stderr, "Too many PE for this problem size => box[0] = -1\n");
141 box[0] = -1;
142 box[1] = -1;
143 box[2] = -1;
144 box[3] = -1;
145 }
146}
147
148#else
149void
150CalcSubSurface(int xmin, int xmax,
151 int ymin, int ymax, int pmin, int pmax, int level, int box[MAXBOX_D], int mype, int pass) {
152 int nbpe = (pmax - pmin + 1);
153 int ny = (int) sqrt(nbpe);
154 int res = (int) (nbpe - ny * ny) / ny;
155 int nx = ny + res;
156 int pex = mype % nx;
157 int pey = mype / nx;
158 int lgx = (xmax - xmin + 1);
159 int incx = lgx / nx;
160 int lgy = (ymax - ymin + 1);
161 int incy = lgy / ny;
162 static int done = 0;
163
164 if (nx * ny != nbpe) {
165 // the closest shape to a square can't be maintain.
166 // Let's find another alternative
167 int divider = 2;
168 int lastdevider = 1;
169 while (divider < (int) sqrt(nbpe)) {
170 if ((nbpe % divider) == 0) {
171 lastdevider = divider;
172 }
173 divider++;
174 }
175
176 // if (mype == 0) printf("Last divider %d\n", lastdevider);
177
178 if(lastdevider == 1) {
179 if (mype == 0) {
180 fprintf(stderr, "\tERROR: %d can't be devided evenly in x and y\n", nbpe);
181 fprintf(stderr, "\tERROR: closest value is %d\n", nx * ny);
182 fprintf(stderr, "\tERROR: please adapt the number of process\n");
183 }
184#ifdef MPI
185#if FTI==0
186 MPI_Finalize();
187#endif
188#if FTI>0
189 FTI_Finalize();
190 MPI_Finalize();
191#endif
192#endif
193 exit(1);
194 }
195 ny = lastdevider;
196 res = (int) (nbpe - ny * ny) / ny;
197 nx = ny + res;
198 pex = mype % nx;
199 pey = mype / nx;
200 incx = lgx / nx;
201 incy = lgy / ny;
202 }
203
204 if ((incx * nx + xmin) < xmax) incx++;
205 if ((incy * ny + ymin) < ymax) incy++;
206
207 if (mype == 0 && !done) {
208 printf("HydroC: Simple decomposition\n");
209 printf("HydroC: nx=%d ny=%d\n", nx, ny);
210 done=1;
211 }
212
213 box[XMIN_D] = pex * incx + xmin;
214 if (box[XMIN_D] < 0) box[XMIN_D] = 0;
215
216 box[XMAX_D] = (pex + 1) * incx + xmin;
217 if (box[XMAX_D] > xmax) box[XMAX_D] = xmax;
218
219
220 box[YMIN_D] = pey * incy + ymin;
221 if (box[YMIN_D] < 0) box[YMIN_D] = 0;
222
223 box[YMAX_D] = (pey + 1) * incy + ymin;
224 if (box[YMAX_D] > ymax) box[YMAX_D] = ymax;
225
226 box[UP_D] = mype + nx;
227 if (box[UP_D] >= nbpe)
228 box[UP_D] = -1;
229 box[DOWN_D] = mype - nx;
230 if (box[DOWN_D] < 0)
231 box[DOWN_D] = -1;
232 box[LEFT_D] = mype - 1;
233 if (pex == 0)
234 box[LEFT_D] = -1;
235 box[RIGHT_D] = mype + 1;
236 if (pex + 1 >= nx)
237 box[RIGHT_D] = -1;
238}
239#endif
240
241#ifdef TESTPGM
242int
243main(int argc, char **argv) {
244 int nx = 16;
245 int ny = 16;
246 int nbproc = 16;
247 int mype = 0;
248 int bord = 0;
249
250 for (mype = 0; mype < nbproc; mype++) {
251 int box[MAXBOX_D] = { -1, -1, -1, -1, -1, -1, -1, -1 };
252 // first pass determin our box
253 CalcSubSurface(0, nx - 1, 0, ny - 1, 0, nbproc - 1, 0, box, mype, 0);
254 // second pass determin our neighbours
255 CalcSubSurface(0, nx - 1, 0, ny - 1, 0, nbproc - 1, 0, box, mype, 1);
256 // printf("[%4d] %4d %4d %4d %4d / %4d %4d %4d %4d \n", mype, box[XMIN_D], box[XMAX_D], box[YMIN_D], box[YMAX_D], box[UP_D], box[DOWN_D], box[LEFT_D], box[RIGHT_D]);
257 if ((box[UP_D] == -1) || (box[DOWN_D] == -1) || (box[LEFT_D] == -1) || (box[RIGHT_D] == -1)
258 )
259 bord++;
260 }
261 printf("B=%d / %d\n", bord, nbproc);
262 return 0;
263}
264#endif
265//EOF
Note: See TracBrowser for help on using the repository browser.