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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 799 - (show annotations)
Mon Aug 7 23:30:53 2006 UTC (13 years, 1 month ago) by ksteube
File MIME type: text/plain
File size: 17173 byte(s)
Moved a #pragma to modify the appropriate line

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26