1 |
|
2 |
/******************************************************* |
3 |
* |
4 |
* Copyright (c) 2003-2008 by University of Queensland |
5 |
* Earth Systems Science Computational Center (ESSCC) |
6 |
* http://www.uq.edu.au/esscc |
7 |
* |
8 |
* Primary Business: Queensland, Australia |
9 |
* Licensed under the Open Software License version 3.0 |
10 |
* http://www.opensource.org/licenses/osl-3.0.php |
11 |
* |
12 |
*******************************************************/ |
13 |
|
14 |
|
15 |
/**************************************************************/ |
16 |
|
17 |
/* Some utility routines: */ |
18 |
|
19 |
/**************************************************************/ |
20 |
|
21 |
#include "Finley.h" |
22 |
#include "Util.h" |
23 |
|
24 |
#ifdef _OPENMP |
25 |
#include <omp.h> |
26 |
#endif |
27 |
|
28 |
/**************************************************************/ |
29 |
|
30 |
/* returns true if any of the values in the short array values is not equalt to Zero */ |
31 |
|
32 |
bool_t Finley_Util_anyNonZeroDouble(dim_t N, double* values) { |
33 |
dim_t q; |
34 |
for (q=0;q<N;++q) if (ABS(values[q])>0) return TRUE; |
35 |
return FALSE; |
36 |
} |
37 |
/**************************************************************/ |
38 |
|
39 |
/* gathers double values out from in by index: */ |
40 |
|
41 |
/* out(1:numData,1:len)=in(1:numData,index(1:len)) */ |
42 |
|
43 |
void Finley_Util_Gather_double(dim_t len,index_t* index,dim_t numData,double* in, double * out){ |
44 |
dim_t s,i; |
45 |
for (s=0;s<len;s++) { |
46 |
for (i=0;i<numData;i++) { |
47 |
out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)]; |
48 |
} |
49 |
} |
50 |
} |
51 |
|
52 |
/**************************************************************/ |
53 |
|
54 |
|
55 |
/* gathers maybelong values out from in by index: */ |
56 |
|
57 |
/* out(1:numData,1:len)=in(1:numData,index(1:len)) */ |
58 |
|
59 |
void Finley_Util_Gather_int(dim_t len,index_t* index,dim_t numData, index_t* in, index_t * out){ |
60 |
dim_t s,i; |
61 |
for (s=0;s<len;s++) { |
62 |
for (i=0;i<numData;i++) { |
63 |
out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)]; |
64 |
} |
65 |
} |
66 |
} |
67 |
|
68 |
/**************************************************************/ |
69 |
|
70 |
/* adds a vector in into out using and index. */ |
71 |
|
72 |
/* out(1:numData,index[p])+=in(1:numData,p) where p = {k=1...len , index[k]<upperBound}*/ |
73 |
|
74 |
|
75 |
void Finley_Util_AddScatter(dim_t len,index_t* index,dim_t numData,double* in,double * out, index_t upperBound){ |
76 |
dim_t i,s; |
77 |
for (s=0;s<len;s++) { |
78 |
for(i=0;i<numData;i++) { |
79 |
if( index[s]<upperBound ) { |
80 |
out[INDEX2(i,index[s],numData)]+=in[INDEX2(i,s,numData)]; |
81 |
} |
82 |
} |
83 |
} |
84 |
} |
85 |
|
86 |
/* multiplies two matrices */ |
87 |
|
88 |
/* A(1:A1,1:A2)=B(1:A1,1:B2)*C(1:B2,1:A2) */ |
89 |
|
90 |
void Finley_Util_SmallMatMult(dim_t A1,dim_t A2, double* A, dim_t B2, double*B, double* C) { |
91 |
dim_t i,j,s; |
92 |
for (i=0;i<A1*A2;i++) A[i]=0; |
93 |
for (i=0;i<A1;i++) { |
94 |
for (j=0;j<A2;j++) { |
95 |
for (s=0;s<B2;s++) { |
96 |
A[INDEX2(i,j,A1)]+=B[INDEX2(i,s,A1)]*C[INDEX2(s,j,B2)]; |
97 |
} |
98 |
} |
99 |
} |
100 |
} |
101 |
|
102 |
/* multiplies a two sets of matries: */ |
103 |
|
104 |
/* A(1:A1,1:A2,i)=B(1:A1,1:B2,i)*C(1:B2,1:A2,i) i=1,len */ |
105 |
|
106 |
void Finley_Util_SmallMatSetMult(dim_t len,dim_t A1,dim_t A2, double* A, dim_t B2, double*B, double* C) { |
107 |
dim_t q,i,j,s; |
108 |
for (i=0;i<A1*A2*len;i++) A[i]=0; |
109 |
for (q=0;q<len;q++) { |
110 |
for (i=0;i<A1;i++) { |
111 |
for (j=0;j<A2;j++) { |
112 |
for (s=0;s<B2;s++) { |
113 |
A[INDEX3(i,j,q,A1,A2)]+=B[INDEX3(i,s,q,A1,B2)]*C[INDEX3(s,j,q,B2,A2)]; |
114 |
} |
115 |
} |
116 |
} |
117 |
} |
118 |
} |
119 |
/* inverts the set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */ |
120 |
/* the determinante is returned. */ |
121 |
|
122 |
void Finley_Util_InvertSmallMat(dim_t len,dim_t dim,double* A,double *invA, double* det){ |
123 |
dim_t q; |
124 |
register double D,A11,A12,A13,A21,A22,A23,A31,A32,A33; |
125 |
|
126 |
switch(dim) { |
127 |
case 1: |
128 |
for (q=0;q<len;q++) { |
129 |
D=A[q]; |
130 |
if (ABS(D) > 0 ){ |
131 |
det[q]=D; |
132 |
D=1./D; |
133 |
invA[q]=D; |
134 |
} else { |
135 |
Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix"); |
136 |
return; |
137 |
} |
138 |
} |
139 |
break; |
140 |
|
141 |
case 2: |
142 |
for (q=0;q<len;q++) { |
143 |
A11=A[INDEX3(0,0,q,2,2)]; |
144 |
A12=A[INDEX3(0,1,q,2,2)]; |
145 |
A21=A[INDEX3(1,0,q,2,2)]; |
146 |
A22=A[INDEX3(1,1,q,2,2)]; |
147 |
|
148 |
D = A11*A22-A12*A21; |
149 |
if (ABS(D) > 0 ){ |
150 |
det[q]=D; |
151 |
D=1./D; |
152 |
invA[INDEX3(0,0,q,2,2)]= A22*D; |
153 |
invA[INDEX3(1,0,q,2,2)]=-A21*D; |
154 |
invA[INDEX3(0,1,q,2,2)]=-A12*D; |
155 |
invA[INDEX3(1,1,q,2,2)]= A11*D; |
156 |
} else { |
157 |
Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix"); |
158 |
return; |
159 |
} |
160 |
} |
161 |
break; |
162 |
|
163 |
case 3: |
164 |
for (q=0;q<len;q++) { |
165 |
A11=A[INDEX3(0,0,q,3,3)]; |
166 |
A21=A[INDEX3(1,0,q,3,3)]; |
167 |
A31=A[INDEX3(2,0,q,3,3)]; |
168 |
A12=A[INDEX3(0,1,q,3,3)]; |
169 |
A22=A[INDEX3(1,1,q,3,3)]; |
170 |
A32=A[INDEX3(2,1,q,3,3)]; |
171 |
A13=A[INDEX3(0,2,q,3,3)]; |
172 |
A23=A[INDEX3(1,2,q,3,3)]; |
173 |
A33=A[INDEX3(2,2,q,3,3)]; |
174 |
|
175 |
D = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22); |
176 |
if (ABS(D) > 0 ){ |
177 |
det[q] =D; |
178 |
D=1./D; |
179 |
invA[INDEX3(0,0,q,3,3)]=(A22*A33-A23*A32)*D; |
180 |
invA[INDEX3(1,0,q,3,3)]=(A31*A23-A21*A33)*D; |
181 |
invA[INDEX3(2,0,q,3,3)]=(A21*A32-A31*A22)*D; |
182 |
invA[INDEX3(0,1,q,3,3)]=(A13*A32-A12*A33)*D; |
183 |
invA[INDEX3(1,1,q,3,3)]=(A11*A33-A31*A13)*D; |
184 |
invA[INDEX3(2,1,q,3,3)]=(A12*A31-A11*A32)*D; |
185 |
invA[INDEX3(0,2,q,3,3)]=(A12*A23-A13*A22)*D; |
186 |
invA[INDEX3(1,2,q,3,3)]=(A13*A21-A11*A23)*D; |
187 |
invA[INDEX3(2,2,q,3,3)]=(A11*A22-A12*A21)*D; |
188 |
} else { |
189 |
Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix"); |
190 |
return; |
191 |
} |
192 |
} |
193 |
break; |
194 |
|
195 |
} |
196 |
return; |
197 |
} |
198 |
|
199 |
/* sets the derterminate of a set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */ |
200 |
|
201 |
void Finley_Util_DetOfSmallMat(dim_t len,dim_t dim,double* A, double* det){ |
202 |
dim_t q; |
203 |
register double A11,A12,A13,A21,A22,A23,A31,A32,A33; |
204 |
|
205 |
switch(dim) { |
206 |
case 1: |
207 |
for (q=0;q<len;q++) { |
208 |
det[q]=A[q]; |
209 |
} |
210 |
break; |
211 |
|
212 |
case 2: |
213 |
for (q=0;q<len;q++) { |
214 |
A11=A[INDEX3(0,0,q,2,2)]; |
215 |
A12=A[INDEX3(0,1,q,2,2)]; |
216 |
A21=A[INDEX3(1,0,q,2,2)]; |
217 |
A22=A[INDEX3(1,1,q,2,2)]; |
218 |
|
219 |
det[q] = A11*A22-A12*A21; |
220 |
} |
221 |
break; |
222 |
|
223 |
case 3: |
224 |
for (q=0;q<len;q++) { |
225 |
A11=A[INDEX3(0,0,q,3,3)]; |
226 |
A21=A[INDEX3(1,0,q,3,3)]; |
227 |
A31=A[INDEX3(2,0,q,3,3)]; |
228 |
A12=A[INDEX3(0,1,q,3,3)]; |
229 |
A22=A[INDEX3(1,1,q,3,3)]; |
230 |
A32=A[INDEX3(2,1,q,3,3)]; |
231 |
A13=A[INDEX3(0,2,q,3,3)]; |
232 |
A23=A[INDEX3(1,2,q,3,3)]; |
233 |
A33=A[INDEX3(2,2,q,3,3)]; |
234 |
|
235 |
det[q] = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22); |
236 |
} |
237 |
break; |
238 |
|
239 |
} |
240 |
return; |
241 |
} |
242 |
/* returns the normalized vector Normal[dim,len] orthogonal to A(:,0,q) and A(:,1,q) in the case of dim=3 */ |
243 |
/* or the vector A(:,0,q) in the case of dim=2 */ |
244 |
|
245 |
void Finley_NormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* Normal) { |
246 |
dim_t q; |
247 |
register double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33,length,invlength; |
248 |
|
249 |
switch(dim) { |
250 |
case 1: |
251 |
for (q=0;q<len;q++) Normal[q] =1; |
252 |
break; |
253 |
case 2: |
254 |
for (q=0;q<len;q++) { |
255 |
A11=A[INDEX3(0,0,q,2,dim1)]; |
256 |
A21=A[INDEX3(1,0,q,2,dim1)]; |
257 |
length = sqrt(A11*A11+A21*A21); |
258 |
if (! length>0) { |
259 |
Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero."); |
260 |
return; |
261 |
} else { |
262 |
invlength=1./length; |
263 |
Normal[INDEX2(0,q,2)]=A21*invlength; |
264 |
Normal[INDEX2(1,q,2)]=-A11*invlength; |
265 |
} |
266 |
} |
267 |
break; |
268 |
case 3: |
269 |
for (q=0;q<len;q++) { |
270 |
A11=A[INDEX3(0,0,q,3,dim1)]; |
271 |
A21=A[INDEX3(1,0,q,3,dim1)]; |
272 |
A31=A[INDEX3(2,0,q,3,dim1)]; |
273 |
A12=A[INDEX3(0,1,q,3,dim1)]; |
274 |
A22=A[INDEX3(1,1,q,3,dim1)]; |
275 |
A32=A[INDEX3(2,1,q,3,dim1)]; |
276 |
CO_A13=A21*A32-A31*A22; |
277 |
CO_A23=A31*A12-A11*A32; |
278 |
CO_A33=A11*A22-A21*A12; |
279 |
length=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33); |
280 |
if (! length>0) { |
281 |
Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero."); |
282 |
return; |
283 |
} else { |
284 |
invlength=1./length; |
285 |
Normal[INDEX2(0,q,3)]=CO_A13*invlength; |
286 |
Normal[INDEX2(1,q,3)]=CO_A23*invlength; |
287 |
Normal[INDEX2(2,q,3)]=CO_A33*invlength; |
288 |
} |
289 |
|
290 |
} |
291 |
break; |
292 |
|
293 |
} |
294 |
return; |
295 |
} |
296 |
|
297 |
/* return the length of the vector which is orthogonal to the vectors A(:,0,q) and A(:,1,q) in the case of dim=3 */ |
298 |
/* or the vector A(:,0,q) in the case of dim=2 */ |
299 |
|
300 |
void Finley_LengthOfNormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* length) { |
301 |
dim_t q; |
302 |
double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33; |
303 |
|
304 |
switch(dim) { |
305 |
case 1: |
306 |
for (q=0;q<len;q++) length[q] =1; |
307 |
break; |
308 |
case 2: |
309 |
for (q=0;q<len;q++) { |
310 |
A11=A[INDEX3(0,0,q,2,dim1)]; |
311 |
A21=A[INDEX3(1,0,q,2,dim1)]; |
312 |
length[q] = sqrt(A11*A11+A21*A21); |
313 |
} |
314 |
break; |
315 |
case 3: |
316 |
for (q=0;q<len;q++) { |
317 |
A11=A[INDEX3(0,0,q,3,dim1)]; |
318 |
A21=A[INDEX3(1,0,q,3,dim1)]; |
319 |
A31=A[INDEX3(2,0,q,3,dim1)]; |
320 |
A12=A[INDEX3(0,1,q,3,dim1)]; |
321 |
A22=A[INDEX3(1,1,q,3,dim1)]; |
322 |
A32=A[INDEX3(2,1,q,3,dim1)]; |
323 |
CO_A13=A21*A32-A31*A22; |
324 |
CO_A23=A31*A12-A11*A32; |
325 |
CO_A33=A11*A22-A21*A12; |
326 |
length[q]=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33); |
327 |
} |
328 |
break; |
329 |
|
330 |
} |
331 |
return; |
332 |
} |
333 |
|
334 |
/* inverts the map map of length len */ |
335 |
/* there is no range checking! */ |
336 |
/* at output Map[invMap[i]]=i for i=0:lenInvMap */ |
337 |
|
338 |
void Finley_Util_InvertMap(dim_t lenInvMap, index_t* invMap,dim_t lenMap, index_t* Map) { |
339 |
dim_t i; |
340 |
for (i=0;i<lenInvMap;i++) invMap[i]=0; |
341 |
for (i=0;i<lenMap;i++) { |
342 |
if (Map[i]>=0) invMap[Map[i]]=i; |
343 |
} |
344 |
} |
345 |
|
346 |
/* orders a Finley_Util_ValueAndIndex array by value */ |
347 |
/* it is assumed that n is large */ |
348 |
|
349 |
int Finley_Util_ValueAndIndex_compar(const void *arg1 , const void *arg2 ) { |
350 |
Finley_Util_ValueAndIndex *e1,*e2; |
351 |
e1=(Finley_Util_ValueAndIndex*) arg1; |
352 |
e2=(Finley_Util_ValueAndIndex*) arg2; |
353 |
if (e1->value < e2->value) return -1; |
354 |
if (e1->value > e2->value) return 1; |
355 |
if (e1->index < e2->index) return -1; |
356 |
if (e1->index > e2->index) return 1; |
357 |
return 0; |
358 |
} |
359 |
|
360 |
void Finley_Util_sortValueAndIndex(dim_t n,Finley_Util_ValueAndIndex* array) { |
361 |
/* OMP : needs parallelization !*/ |
362 |
qsort(array,n,sizeof(Finley_Util_ValueAndIndex),Finley_Util_ValueAndIndex_compar); |
363 |
} |
364 |
|
365 |
|
366 |
/**************************************************************/ |
367 |
|
368 |
/* calculates the minimum value from a dim X N integer array */ |
369 |
|
370 |
index_t Finley_Util_getMinInt(dim_t dim,dim_t N,index_t* values) { |
371 |
dim_t i,j; |
372 |
index_t out,out_local; |
373 |
out=INDEX_T_MAX; |
374 |
if (values!=NULL && dim*N>0 ) { |
375 |
out=values[0]; |
376 |
#pragma omp parallel private(out_local) |
377 |
{ |
378 |
out_local=out; |
379 |
#pragma omp for private(i,j) schedule(static) |
380 |
for (j=0;j<N;j++) { |
381 |
for (i=0;i<dim;i++) out_local=MIN(out_local,values[INDEX2(i,j,dim)]); |
382 |
} |
383 |
#pragma omp critical |
384 |
out=MIN(out_local,out); |
385 |
} |
386 |
} |
387 |
return out; |
388 |
} |
389 |
|
390 |
/* calculates the maximum value from a dim X N integer array */ |
391 |
|
392 |
index_t Finley_Util_getMaxInt(dim_t dim,dim_t N,index_t* values) { |
393 |
dim_t i,j; |
394 |
index_t out,out_local; |
395 |
out=-INDEX_T_MAX; |
396 |
if (values!=NULL && dim*N>0 ) { |
397 |
out=values[0]; |
398 |
#pragma omp parallel private(out_local) |
399 |
{ |
400 |
out_local=out; |
401 |
#pragma omp for private(i,j) schedule(static) |
402 |
for (j=0;j<N;j++) { |
403 |
for (i=0;i<dim;i++) out_local=MAX(out_local,values[INDEX2(i,j,dim)]); |
404 |
} |
405 |
#pragma omp critical |
406 |
out=MAX(out_local,out); |
407 |
} |
408 |
} |
409 |
return out; |
410 |
} |
411 |
/**************************************************************/ |
412 |
|
413 |
/* calculates the minimum value from a dim X N integer array */ |
414 |
|
415 |
index_t Finley_Util_getFlaggedMinInt(dim_t dim,dim_t N,index_t* values, index_t ignore) { |
416 |
dim_t i,j; |
417 |
index_t out,out_local; |
418 |
out=INDEX_T_MAX; |
419 |
if (values!=NULL && dim*N>0 ) { |
420 |
out=values[0]; |
421 |
#pragma omp parallel private(out_local) |
422 |
{ |
423 |
out_local=out; |
424 |
#pragma omp for private(i,j) schedule(static) |
425 |
for (j=0;j<N;j++) { |
426 |
for (i=0;i<dim;i++) if (values[INDEX2(i,j,dim)]!=ignore) out_local=MIN(out_local,values[INDEX2(i,j,dim)]); |
427 |
} |
428 |
#pragma omp critical |
429 |
out=MIN(out_local,out); |
430 |
} |
431 |
} |
432 |
return out; |
433 |
} |
434 |
|
435 |
/* calculates the maximum value from a dim X N integer array */ |
436 |
|
437 |
index_t Finley_Util_getFlaggedMaxInt(dim_t dim,dim_t N,index_t* values, index_t ignore) { |
438 |
dim_t i,j; |
439 |
index_t out,out_local; |
440 |
out=-INDEX_T_MAX; |
441 |
if (values!=NULL && dim*N>0 ) { |
442 |
out=values[0]; |
443 |
#pragma omp parallel private(out_local) |
444 |
{ |
445 |
out_local=out; |
446 |
#pragma omp for private(i,j) schedule(static) |
447 |
for (j=0;j<N;j++) { |
448 |
for (i=0;i<dim;i++) if (values[INDEX2(i,j,dim)]!=ignore) out_local=MAX(out_local,values[INDEX2(i,j,dim)]); |
449 |
} |
450 |
#pragma omp critical |
451 |
out=MAX(out_local,out); |
452 |
} |
453 |
} |
454 |
return out; |
455 |
} |
456 |
|
457 |
/* set the index of the positive entries in mask. The length of index is returned. */ |
458 |
|
459 |
dim_t Finley_Util_packMask(dim_t N,index_t* mask,index_t* index) { |
460 |
dim_t out,k; |
461 |
out=0; |
462 |
/*OMP */ |
463 |
for (k=0;k<N;k++) { |
464 |
if (mask[k]>=0) { |
465 |
index[out]=k; |
466 |
out++; |
467 |
} |
468 |
} |
469 |
return out; |
470 |
} |
471 |
|
472 |
/* returns true if array contains value */ |
473 |
bool_t Finley_Util_isAny(dim_t N,index_t* array,index_t value) { |
474 |
bool_t out=FALSE; |
475 |
dim_t i; |
476 |
#pragma omp parallel for private(i) schedule(static) reduction(||:out) |
477 |
for (i=0;i<N;i++) out = out || (array[i]==value); |
478 |
return out; |
479 |
} |
480 |
/* calculates the cummultative sum in array and returns the total sum */ |
481 |
index_t Finley_Util_cumsum(dim_t N,index_t* array) { |
482 |
index_t out=0,tmp; |
483 |
dim_t i; |
484 |
#ifdef _OPENMP |
485 |
index_t *partial_sums=NULL, sum; |
486 |
partial_sums=TMPMEMALLOC(omp_get_max_threads(),index_t); |
487 |
#pragma omp parallel private(sum,i,tmp) |
488 |
{ |
489 |
sum=0; |
490 |
#pragma omp for schedule(static) |
491 |
for (i=0;i<N;++i) sum+=array[i]; |
492 |
partial_sums[omp_get_thread_num()]=sum; |
493 |
#pragma omp barrier |
494 |
#pragma omp master |
495 |
{ |
496 |
out=0; |
497 |
for (i=0;i<omp_get_max_threads();++i) { |
498 |
tmp=out; |
499 |
out+=partial_sums[i]; |
500 |
partial_sums[i]=tmp; |
501 |
} |
502 |
} |
503 |
#pragma omp barrier |
504 |
sum=partial_sums[omp_get_thread_num()]; |
505 |
#pragma omp for schedule(static) |
506 |
for (i=0;i<N;++i) { |
507 |
tmp=sum; |
508 |
sum+=array[i]; |
509 |
array[i]=tmp; |
510 |
} |
511 |
} |
512 |
TMPMEMFREE(partial_sums); |
513 |
#else |
514 |
for (i=0;i<N;++i) { |
515 |
tmp=out; |
516 |
out+=array[i]; |
517 |
array[i]=tmp; |
518 |
} |
519 |
#endif |
520 |
return out; |
521 |
} |
522 |
void Finley_Util_setValuesInUse(const index_t *values, const dim_t numValues, dim_t *numValuesInUse, index_t **valuesInUse, Paso_MPIInfo* mpiinfo) |
523 |
{ |
524 |
dim_t i; |
525 |
index_t lastFoundValue=INDEX_T_MIN, minFoundValue, local_minFoundValue, *newValuesInUse=NULL; |
526 |
register index_t itmp; |
527 |
bool_t allFound=FALSE; |
528 |
dim_t nv=0; |
529 |
|
530 |
while (! allFound) { |
531 |
/* |
532 |
* find smallest value bigger than lastFoundValue |
533 |
*/ |
534 |
minFoundValue=INDEX_T_MAX; |
535 |
#pragma omp parallel private(local_minFoundValue) |
536 |
{ |
537 |
local_minFoundValue=minFoundValue; |
538 |
#pragma omp for private(i,itmp) schedule(static) |
539 |
for (i=0;i< numValues;i++) { |
540 |
itmp=values[i]; |
541 |
if ((itmp>lastFoundValue) && (itmp<local_minFoundValue)) local_minFoundValue=itmp; |
542 |
} |
543 |
#pragma omp critical |
544 |
minFoundValue=MIN(local_minFoundValue,minFoundValue); |
545 |
} |
546 |
#ifdef PASO_MPI |
547 |
local_minFoundValue=minFoundValue; |
548 |
MPI_Allreduce(&local_minFoundValue,&minFoundValue, 1, MPI_INT, MPI_MAX, mpiinfo->comm ); |
549 |
#endif |
550 |
|
551 |
/* if we found a new tag we need to add this too the valuesInUseList */ |
552 |
|
553 |
if (minFoundValue < INDEX_T_MAX) { |
554 |
newValuesInUse=MEMALLOC(nv+1,index_t); |
555 |
if (*valuesInUse!=NULL) { |
556 |
memcpy(newValuesInUse,*valuesInUse,sizeof(index_t)*nv); |
557 |
MEMFREE(*valuesInUse); |
558 |
} |
559 |
newValuesInUse[nv]=minFoundValue; |
560 |
*valuesInUse=newValuesInUse; |
561 |
newValuesInUse=NULL; |
562 |
nv++; |
563 |
lastFoundValue=minFoundValue; |
564 |
} else { |
565 |
allFound=TRUE; |
566 |
} |
567 |
} |
568 |
*numValuesInUse=nv; |
569 |
} |
570 |
|
571 |
|
572 |
#ifdef PASO_MPI |
573 |
void Finley_printDoubleArray( FILE *fid, dim_t n, double *array, char *name ) |
574 |
{ |
575 |
index_t i; |
576 |
|
577 |
if( name ) |
578 |
fprintf( fid, "%s [ ", name ); |
579 |
else |
580 |
fprintf( fid, "[ " ); |
581 |
for( i=0; i<(n<60 ? n : 60); i++ ) |
582 |
fprintf( fid, "%g ", array[i] ); |
583 |
if( n>=30 ) |
584 |
fprintf( fid, "... " ); |
585 |
fprintf( fid, "]\n" ); |
586 |
} |
587 |
void Finley_printIntArray( FILE *fid, dim_t n, int *array, char *name ) |
588 |
{ |
589 |
index_t i; |
590 |
|
591 |
if( name ) |
592 |
fprintf( fid, "%s [ ", name ); |
593 |
else |
594 |
fprintf( fid, "[ " ); |
595 |
for( i=0; i<(n<60 ? n : 60); i++ ) |
596 |
fprintf( fid, "%d ", array[i] ); |
597 |
if( n>=30 ) |
598 |
fprintf( fid, "... " ); |
599 |
fprintf( fid, "]\n" ); |
600 |
} |
601 |
void Finley_printMaskArray( FILE *fid, dim_t n, int *array, char *name ) |
602 |
{ |
603 |
index_t i; |
604 |
|
605 |
if( name ) |
606 |
fprintf( fid, "%s [ ", name ); |
607 |
else |
608 |
fprintf( fid, "[ " ); |
609 |
for( i=0; i<(n<60 ? n : 60); i++ ) |
610 |
if( array[i]!=-1 ) |
611 |
fprintf( fid, "%3d ", array[i] ); |
612 |
else |
613 |
fprintf( fid, " * " ); |
614 |
if( n>=30 ) |
615 |
fprintf( fid, "... " ); |
616 |
fprintf( fid, "]\n" ); |
617 |
} |
618 |
#endif |