/[escript]/trunk/finley/src/Util.c
ViewVC logotype

Contents of /trunk/finley/src/Util.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 751 - (show annotations)
Mon Jun 26 01:46:34 2006 UTC (12 years, 10 months ago) by bcumming
File MIME type: text/plain
File size: 17647 byte(s)
Changes relating to the MPI version of escript
The standard OpenMP version of escript is unchanged

- updated data types (Finley_Mesh, Finley_NodeFile, etc) to store meshes
  over multiple MPI processes.
- added CommBuffer code in Paso for communication of Data associated
  with distributed meshes
- updates in Finley and Escript to support distributed data and operations
  on distributed data (such as interpolation).
- construction of RHS in MPI, so that simple explicit schemes (such as
  /docs/examples/wave.py without IO and the Locator) can run in MPI.
- updated mesh generation for first order line, rectangle and brick
  meshes and second order line meshes in MPI.        
- small changes to trunk/SConstruct and trunk/scons/ess_options.py to
  build the MPI version, these changes are turned off by default.

1 /*
2 ************************************************************
3 * Copyright 2006 by ACcESS MNRF *
4 * *
5 * http://www.access.edu.au *
6 * Primary Business: Queensland, Australia *
7 * Licensed under the Open Software License version 3.0 *
8 * http://www.opensource.org/licenses/osl-3.0.php *
9 * *
10 ************************************************************
11 */
12
13 /**************************************************************/
14
15 /* Some utility routines: */
16
17 /**************************************************************/
18
19 /* author: gross@access.edu.au */
20 /* Version: $Id$ */
21
22 /**************************************************************/
23
24 #include "Finley.h"
25 #include "Util.h"
26 #ifdef _OPENMP
27 #include <omp.h>
28 #endif
29
30 /**************************************************************/
31
32 /* returns true if any of the values in the short array values is not equalt to Zero */
33
34 bool_t Finley_Util_anyNonZeroDouble(dim_t N, double* values) {
35 dim_t q;
36 for (q=0;q<N;++q) if (ABS(values[q])>0) return TRUE;
37 return FALSE;
38 }
39 /**************************************************************/
40
41 /* gathers double values out from in by index: */
42
43 /* out(1:numData,1:len)=in(1:numData,index(1:len)) */
44
45 void Finley_Util_Gather_double(dim_t len,index_t* index,dim_t numData,double* in, double * out){
46 dim_t s,i;
47 for (s=0;s<len;s++) {
48 for (i=0;i<numData;i++) {
49 out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)];
50 }
51 }
52 }
53
54 /**************************************************************/
55
56
57 /* gathers maybelong values out from in by index: */
58
59 /* out(1:numData,1:len)=in(1:numData,index(1:len)) */
60
61 void Finley_Util_Gather_int(dim_t len,index_t* index,dim_t numData, index_t* in, index_t * out){
62 dim_t s,i;
63 for (s=0;s<len;s++) {
64 for (i=0;i<numData;i++) {
65 out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)];
66 }
67 }
68 }
69
70 /**************************************************************/
71
72 /* adds a vector in into out using and index. */
73
74 /* out(1:numData,index(1:len))+=in(1:numData,1:len) */
75
76 void Finley_Util_AddScatter(dim_t len,index_t* index,dim_t numData,double* in,double * out){
77 dim_t i,s;
78 for (s=0;s<len;s++) {
79 for(i=0;i<numData;i++) {
80 #pragma omp atomic
81 out[INDEX2(i,index[s],numData)]+=in[INDEX2(i,s,numData)];
82 }
83 }
84 }
85
86 #ifdef PASO_MPI
87 /* same as AddScatter(), but checks that value index[] is below an upper bound upperBound before
88 addition. This is used to ensure that only the influence of local DOF is added */
89 /* out(1:numData,index[p])+=in(1:numData,p)
90 where p = {k=1...len , index[k]<upperBound}*/
91 void Finley_Util_AddScatter_upperBound(dim_t len,index_t* index,dim_t numData,double* in,double * out, index_t upperBound){
92 dim_t i,s;
93 for (s=0;s<len;s++) {
94 for(i=0;i<numData;i++) {
95 //#pragma omp atomic
96 if( index[s]<upperBound )
97 out[INDEX2(i,index[s],numData)]+=in[INDEX2(i,s,numData)];
98 }
99 }
100 }
101
102
103 #endif
104
105 /* multiplies two matrices */
106
107 /* A(1:A1,1:A2)=B(1:A1,1:B2)*C(1:B2,1:A2) */
108
109 void Finley_Util_SmallMatMult(dim_t A1,dim_t A2, double* A, dim_t B2, double*B, double* C) {
110 dim_t i,j,s;
111 for (i=0;i<A1*A2;i++) A[i]=0;
112 for (i=0;i<A1;i++) {
113 for (j=0;j<A2;j++) {
114 for (s=0;s<B2;s++) {
115 A[INDEX2(i,j,A1)]+=B[INDEX2(i,s,A1)]*C[INDEX2(s,j,B2)];
116 }
117 }
118 }
119 }
120
121 /* multiplies a two sets of matries: */
122
123 /* A(1:A1,1:A2,i)=B(1:A1,1:B2,i)*C(1:B2,1:A2,i) i=1,len */
124
125 void Finley_Util_SmallMatSetMult(dim_t len,dim_t A1,dim_t A2, double* A, dim_t B2, double*B, double* C) {
126 dim_t q,i,j,s;
127 for (i=0;i<A1*A2*len;i++) A[i]=0;
128 for (q=0;q<len;q++) {
129 for (i=0;i<A1;i++) {
130 for (j=0;j<A2;j++) {
131 for (s=0;s<B2;s++) {
132 A[INDEX3(i,j,q,A1,A2)]+=B[INDEX3(i,s,q,A1,B2)]*C[INDEX3(s,j,q,B2,A2)];
133 }
134 }
135 }
136 }
137 }
138 /* inverts the set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */
139 /* the determinante is returned. */
140
141 void Finley_Util_InvertSmallMat(dim_t len,dim_t dim,double* A,double *invA, double* det){
142 dim_t q;
143 register double D,A11,A12,A13,A21,A22,A23,A31,A32,A33;
144
145 switch(dim) {
146 case 1:
147 for (q=0;q<len;q++) {
148 D=A[q];
149 if (ABS(D) > 0 ){
150 det[q]=D;
151 D=1./D;
152 invA[q]=D;
153 } else {
154 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
155 return;
156 }
157 }
158 break;
159
160 case 2:
161 for (q=0;q<len;q++) {
162 A11=A[INDEX3(0,0,q,2,2)];
163 A12=A[INDEX3(0,1,q,2,2)];
164 A21=A[INDEX3(1,0,q,2,2)];
165 A22=A[INDEX3(1,1,q,2,2)];
166
167 D = A11*A22-A12*A21;
168 if (ABS(D) > 0 ){
169 det[q]=D;
170 D=1./D;
171 invA[INDEX3(0,0,q,2,2)]= A22*D;
172 invA[INDEX3(1,0,q,2,2)]=-A21*D;
173 invA[INDEX3(0,1,q,2,2)]=-A12*D;
174 invA[INDEX3(1,1,q,2,2)]= A11*D;
175 } else {
176 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
177 return;
178 }
179 }
180 break;
181
182 case 3:
183 for (q=0;q<len;q++) {
184 A11=A[INDEX3(0,0,q,3,3)];
185 A21=A[INDEX3(1,0,q,3,3)];
186 A31=A[INDEX3(2,0,q,3,3)];
187 A12=A[INDEX3(0,1,q,3,3)];
188 A22=A[INDEX3(1,1,q,3,3)];
189 A32=A[INDEX3(2,1,q,3,3)];
190 A13=A[INDEX3(0,2,q,3,3)];
191 A23=A[INDEX3(1,2,q,3,3)];
192 A33=A[INDEX3(2,2,q,3,3)];
193
194 D = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22);
195 if (ABS(D) > 0 ){
196 det[q] =D;
197 D=1./D;
198 invA[INDEX3(0,0,q,3,3)]=(A22*A33-A23*A32)*D;
199 invA[INDEX3(1,0,q,3,3)]=(A31*A23-A21*A33)*D;
200 invA[INDEX3(2,0,q,3,3)]=(A21*A32-A31*A22)*D;
201 invA[INDEX3(0,1,q,3,3)]=(A13*A32-A12*A33)*D;
202 invA[INDEX3(1,1,q,3,3)]=(A11*A33-A31*A13)*D;
203 invA[INDEX3(2,1,q,3,3)]=(A12*A31-A11*A32)*D;
204 invA[INDEX3(0,2,q,3,3)]=(A12*A23-A13*A22)*D;
205 invA[INDEX3(1,2,q,3,3)]=(A13*A21-A11*A23)*D;
206 invA[INDEX3(2,2,q,3,3)]=(A11*A22-A12*A21)*D;
207 } else {
208 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
209 return;
210 }
211 }
212 break;
213
214 }
215 return;
216 }
217
218 /* sets the derterminate of a set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */
219
220 void Finley_Util_DetOfSmallMat(dim_t len,dim_t dim,double* A, double* det){
221 dim_t q;
222 register double A11,A12,A13,A21,A22,A23,A31,A32,A33;
223
224 switch(dim) {
225 case 1:
226 for (q=0;q<len;q++) {
227 det[q]=A[q];
228 }
229 break;
230
231 case 2:
232 for (q=0;q<len;q++) {
233 A11=A[INDEX3(0,0,q,2,2)];
234 A12=A[INDEX3(0,1,q,2,2)];
235 A21=A[INDEX3(1,0,q,2,2)];
236 A22=A[INDEX3(1,1,q,2,2)];
237
238 det[q] = A11*A22-A12*A21;
239 }
240 break;
241
242 case 3:
243 for (q=0;q<len;q++) {
244 A11=A[INDEX3(0,0,q,3,3)];
245 A21=A[INDEX3(1,0,q,3,3)];
246 A31=A[INDEX3(2,0,q,3,3)];
247 A12=A[INDEX3(0,1,q,3,3)];
248 A22=A[INDEX3(1,1,q,3,3)];
249 A32=A[INDEX3(2,1,q,3,3)];
250 A13=A[INDEX3(0,2,q,3,3)];
251 A23=A[INDEX3(1,2,q,3,3)];
252 A33=A[INDEX3(2,2,q,3,3)];
253
254 det[q] = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22);
255 }
256 break;
257
258 }
259 return;
260 }
261 /* returns the normalized vector Normal[dim,len] orthogonal to A(:,0,q) and A(:,1,q) in the case of dim=3 */
262 /* or the vector A(:,0,q) in the case of dim=2 */
263
264 void Finley_NormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* Normal) {
265 dim_t q;
266 register double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33,length,invlength;
267
268 switch(dim) {
269 case 1:
270 for (q=0;q<len;q++) Normal[q] =1;
271 break;
272 case 2:
273 for (q=0;q<len;q++) {
274 A11=A[INDEX3(0,0,q,2,dim1)];
275 A21=A[INDEX3(1,0,q,2,dim1)];
276 length = sqrt(A11*A11+A21*A21);
277 if (! length>0) {
278 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero.");
279 return;
280 } else {
281 invlength=1./length;
282 Normal[INDEX2(0,q,2)]=A21*invlength;
283 Normal[INDEX2(1,q,2)]=-A11*invlength;
284 }
285 }
286 break;
287 case 3:
288 for (q=0;q<len;q++) {
289 A11=A[INDEX3(0,0,q,3,dim1)];
290 A21=A[INDEX3(1,0,q,3,dim1)];
291 A31=A[INDEX3(2,0,q,3,dim1)];
292 A12=A[INDEX3(0,1,q,3,dim1)];
293 A22=A[INDEX3(1,1,q,3,dim1)];
294 A32=A[INDEX3(2,1,q,3,dim1)];
295 CO_A13=A21*A32-A31*A22;
296 CO_A23=A31*A12-A11*A32;
297 CO_A33=A11*A22-A21*A12;
298 length=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33);
299 if (! length>0) {
300 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero.");
301 return;
302 } else {
303 invlength=1./length;
304 Normal[INDEX2(0,q,3)]=CO_A13*invlength;
305 Normal[INDEX2(1,q,3)]=CO_A23*invlength;
306 Normal[INDEX2(2,q,3)]=CO_A33*invlength;
307 }
308
309 }
310 break;
311
312 }
313 return;
314 }
315
316 /* 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 */
317 /* or the vector A(:,0,q) in the case of dim=2 */
318
319 void Finley_LengthOfNormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* length) {
320 dim_t q;
321 double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33;
322
323 switch(dim) {
324 case 1:
325 for (q=0;q<len;q++) length[q] =1;
326 break;
327 case 2:
328 for (q=0;q<len;q++) {
329 A11=A[INDEX3(0,0,q,2,dim1)];
330 A21=A[INDEX3(1,0,q,2,dim1)];
331 length[q] = sqrt(A11*A11+A21*A21);
332 }
333 break;
334 case 3:
335 for (q=0;q<len;q++) {
336 A11=A[INDEX3(0,0,q,3,dim1)];
337 A21=A[INDEX3(1,0,q,3,dim1)];
338 A31=A[INDEX3(2,0,q,3,dim1)];
339 A12=A[INDEX3(0,1,q,3,dim1)];
340 A22=A[INDEX3(1,1,q,3,dim1)];
341 A32=A[INDEX3(2,1,q,3,dim1)];
342 CO_A13=A21*A32-A31*A22;
343 CO_A23=A31*A12-A11*A32;
344 CO_A33=A11*A22-A21*A12;
345 length[q]=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33);
346 }
347 break;
348
349 }
350 return;
351 }
352
353 /* inverts the map map of length len */
354 /* there is no range checking! */
355 /* at output Map[invMap[i]]=i for i=0:lenInvMap */
356
357 void Finley_Util_InvertMap(dim_t lenInvMap, index_t* invMap,dim_t lenMap, index_t* Map) {
358 dim_t i;
359 for (i=0;i<lenInvMap;i++) invMap[i]=0;
360 for (i=0;i<lenMap;i++) {
361 if (Map[i]>=0) invMap[Map[i]]=i;
362 }
363 }
364
365 /* orders a Finley_Util_ValueAndIndex array by value */
366 /* it is assumed that n is large */
367
368 int Finley_Util_ValueAndIndex_compar(const void *arg1 , const void *arg2 ) {
369 Finley_Util_ValueAndIndex *e1,*e2;
370 e1=(Finley_Util_ValueAndIndex*) arg1;
371 e2=(Finley_Util_ValueAndIndex*) arg2;
372 if (e1->value < e2->value) return -1;
373 if (e1->value > e2->value) return 1;
374 return 0;
375 }
376 void Finley_Util_sortValueAndIndex(dim_t n,Finley_Util_ValueAndIndex* array) {
377 /* OMP : needs parallelization !*/
378 qsort(array,n,sizeof(Finley_Util_ValueAndIndex),Finley_Util_ValueAndIndex_compar);
379 }
380
381
382 /**************************************************************/
383
384 /* calculates the minimum value from a dim X N integer array */
385
386 index_t Finley_Util_getMinInt(dim_t dim,dim_t N,index_t* values) {
387 dim_t i,j;
388 index_t out,out_local;
389 out=INDEX_T_MAX;
390 if (values!=NULL && dim*N>0 ) {
391 out=values[0];
392 #pragma omp parallel private(out_local)
393 {
394 out_local=out;
395 #pragma omp for private(i,j) schedule(static)
396 for (j=0;j<N;j++) {
397 for (i=0;i<dim;i++) out_local=MIN(out_local,values[INDEX2(i,j,dim)]);
398 }
399 #pragma omp critical
400 out=MIN(out_local,out);
401 }
402 }
403 return out;
404 }
405
406 /* calculates the maximum value from a dim X N integer array */
407
408 index_t Finley_Util_getMaxInt(dim_t dim,dim_t N,index_t* values) {
409 dim_t i,j;
410 index_t out,out_local;
411 out=-INDEX_T_MAX;
412 if (values!=NULL && dim*N>0 ) {
413 out=values[0];
414 #pragma omp parallel private(out_local)
415 {
416 out_local=out;
417 #pragma omp for private(i,j) schedule(static)
418 for (j=0;j<N;j++) {
419 for (i=0;i<dim;i++) out_local=MAX(out_local,values[INDEX2(i,j,dim)]);
420 }
421 #pragma omp critical
422 out=MAX(out_local,out);
423 }
424 }
425 return out;
426 }
427
428 /* set the index of the positive entries in mask. The length of index is returned. */
429
430 dim_t Finley_Util_packMask(dim_t N,index_t* mask,index_t* index) {
431 dim_t out,k;
432 out=0;
433 /*OMP */
434 for (k=0;k<N;k++) {
435 if (mask[k]>=0) {
436 index[out]=k;
437 out++;
438 }
439 }
440 return out;
441 }
442
443 /* returns true if array contains value */
444 bool_t Finley_Util_isAny(dim_t N,index_t* array,index_t value) {
445 bool_t out=FALSE;
446 dim_t i;
447 #pragma omp parallel for private(i) schedule(static) reduction(||:out)
448 for (i=0;i<N;i++) out = out || (array[i]==value);
449 return out;
450 }
451 /* calculates the cummultative sum in array and returns the total sum */
452 index_t Finley_Util_cumsum(dim_t N,index_t* array) {
453 index_t out=0,tmp;
454 dim_t i;
455 #ifdef _OPENMP
456 index_t partial_sums[omp_get_max_threads()],sum;
457 #pragma omp parallel private(sum,i,tmp)
458 {
459 sum=0;
460 #pragma omp for schedule(static)
461 for (i=0;i<N;++i) sum+=array[i];
462 partial_sums[omp_get_thread_num()]=sum;
463 #pragma omp barrier
464 #pragma omp master
465 {
466 out=0;
467 for (i=0;i<omp_get_max_threads();++i) {
468 tmp=out;
469 out+=partial_sums[i];
470 partial_sums[i]=tmp;
471 }
472 }
473 #pragma omp barrier
474 sum=partial_sums[omp_get_thread_num()];
475 #pragma omp for schedule(static)
476 for (i=0;i<N;++i) {
477 tmp=sum;
478 sum+=array[i];
479 array[i]=tmp;
480 }
481 }
482 #else
483 for (i=0;i<N;++i) {
484 tmp=out;
485 out+=array[i];
486 array[i]=tmp;
487 }
488 #endif
489 return out;
490 }
491
492 void Finley_copyDouble(dim_t n,double* source, double* target) {
493 dim_t i;
494 for (i=0;i<n;i++) target[i]=source[i];
495 }
496
497 #ifdef PASO_MPI
498 void Finley_printDoubleArray( FILE *fid, dim_t n, double *array, char *name )
499 {
500 index_t i;
501
502 if( name )
503 fprintf( fid, "%s [ ", name );
504 else
505 fprintf( fid, "[ " );
506 for( i=0; i<(n<30 ? n : 30); i++ )
507 fprintf( fid, "%g ", array[i] );
508 if( n>=30 )
509 fprintf( fid, "... " );
510 fprintf( fid, "]\n" );
511 }
512 void Finley_printIntArray( FILE *fid, dim_t n, int *array, char *name )
513 {
514 index_t i;
515
516 if( name )
517 fprintf( fid, "%s [ ", name );
518 else
519 fprintf( fid, "[ " );
520 for( i=0; i<(n<30 ? n : 30); i++ )
521 fprintf( fid, "%d ", array[i] );
522 if( n>=30 )
523 fprintf( fid, "... " );
524 fprintf( fid, "]\n" );
525 }
526 void Finley_printMaskArray( FILE *fid, dim_t n, int *array, char *name )
527 {
528 index_t i;
529
530 if( name )
531 fprintf( fid, "%s [ ", name );
532 else
533 fprintf( fid, "[ " );
534 for( i=0; i<(n<30 ? n : 30); i++ )
535 if( array[i]!=-1 )
536 fprintf( fid, "%d ", array[i] );
537 else
538 fprintf( fid, "* " );
539 if( n>=30 )
540 fprintf( fid, "... " );
541 fprintf( fid, "]\n" );
542 }
543 #endif
544
545 /*
546 * Revision 1.8 2005/08/12 01:45:43 jgs
547 *
548 * Revision 1.7.2.2 2005/09/07 06:26:22 gross
549 * the solver from finley are put into the standalone package paso now
550 *
551 * Revision 1.7.2.1 2005/08/04 22:41:11 gross
552 * some extra routines for finley that might speed-up RHS assembling in some cases (not actived right now)
553 *
554 * Revision 1.7 2005/07/08 04:07:59 jgs
555 * Merge of development branch back to main trunk on 2005-07-08
556 *
557 * Revision 1.1.1.1.2.4 2005/06/29 02:34:57 gross
558 * some changes towards 64 integers in finley
559 *
560 * Revision 1.1.1.1.2.3 2005/03/02 23:35:06 gross
561 * reimplementation of the ILU in Finley. block size>1 still needs some testing
562 *
563 * Revision 1.1.1.1.2.2 2005/02/18 02:27:31 gross
564 * two function that will be used for a reimplementation of the ILU preconditioner
565 *
566 * Revision 1.1.1.1.2.1 2004/11/12 06:58:19 gross
567 * a lot of changes to get the linearPDE class running: most important change is that there is no matrix format exposed to the user anymore. the format is chosen by the Domain according to the solver and symmetry
568 *
569 * Revision 1.1.1.1 2004/10/26 06:53:57 jgs
570 * initial import of project esys2
571 *
572 * Revision 1.3 2004/08/26 12:03:52 gross
573 * Some other bug in Finley_Assemble_gradient fixed.
574 *
575 * Revision 1.2 2004/07/02 04:21:13 gross
576 * Finley C code has been included
577 *
578 * Revision 1.1.1.1 2004/06/24 04:00:40 johng
579 * Initial version of eys using boost-python.
580 *
581 *
582 */

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26