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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 763 - (show annotations)
Fri Jun 30 05:52:31 2006 UTC (13 years, 3 months ago) by gross
File MIME type: text/plain
File size: 17731 byte(s)
some modification in the element reordering which should solve some problems 
caused by the fact that qsort may handels == items differently on different platform.


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26