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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 616 - (show annotations)
Wed Mar 22 02:46:56 2006 UTC (13 years, 6 months ago) by elspeth
File MIME type: text/plain
File size: 16008 byte(s)
Copyright added to more source files.

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26