/[escript]/branches/domexper/dudley/src/Util.c
ViewVC logotype

Contents of /branches/domexper/dudley/src/Util.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 798 - (show annotations)
Fri Aug 4 01:05:36 2006 UTC (13 years, 6 months ago) by gross
Original Path: trunk/finley/src/Util.c
File MIME type: text/plain
File size: 17164 byte(s)
Reimplementation of the assemblage with persistent jacobeans.
There are also a few changes to the tests which has now
dramatically reduced the memory demand.


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26