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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 471 - (show annotations)
Fri Jan 27 01:33:02 2006 UTC (13 years, 7 months ago) by jgs
File MIME type: text/plain
File size: 16330 byte(s)
reorganise finley src tree to remove inc dir and src/finley directory

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26