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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2748 - (show annotations)
Tue Nov 17 07:32:59 2009 UTC (10 years ago) by gross
File MIME type: text/plain
File size: 19205 byte(s)
Macro elements are implemented now. VTK writer for macro elements still needs testing.
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2009 by University of Queensland
5 * Earth Systems Science Computational Center (ESSCC)
6 * http://www.uq.edu.au/esscc
7 *
8 * Primary Business: Queensland, Australia
9 * Licensed under the Open Software License version 3.0
10 * http://www.opensource.org/licenses/osl-3.0.php
11 *
12 *******************************************************/
13
14
15 /**************************************************************/
16
17 /* Some utility routines: */
18
19 /**************************************************************/
20
21 #include "Finley.h"
22 #include "Util.h"
23
24 #ifdef _OPENMP
25 #include <omp.h>
26 #endif
27
28 /**************************************************************/
29
30 /* returns true if any of the values in the short array values is not equalt to Zero */
31
32 bool_t Finley_Util_anyNonZeroDouble(dim_t N, double* values) {
33 dim_t q;
34 for (q=0;q<N;++q) if (ABS(values[q])>0) return TRUE;
35 return FALSE;
36 }
37 /**************************************************************/
38
39 /* gathers double values out from in by index: */
40
41 /* out(1:numData,1:len)=in(1:numData,index(1:len)) */
42
43 void Finley_Util_Gather_double(dim_t len,index_t* index,dim_t numData,double* in, double * out){
44 dim_t s,i;
45 for (s=0;s<len;s++) {
46 for (i=0;i<numData;i++) {
47 out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)];
48 }
49 }
50 }
51
52 /**************************************************************/
53
54
55 /* gathers maybelong values out from in by index: */
56
57 /* out(1:numData,1:len)=in(1:numData,index(1:len)) */
58
59 void Finley_Util_Gather_int(dim_t len,index_t* index,dim_t numData, index_t* in, index_t * out){
60 dim_t s,i;
61 for (s=0;s<len;s++) {
62 for (i=0;i<numData;i++) {
63 out[INDEX2(i,s,numData)]=in[INDEX2(i,index[s],numData)];
64 }
65 }
66 }
67
68 /**************************************************************/
69
70 /* adds a vector in into out using and index. */
71
72 /* out(1:numData,index[p])+=in(1:numData,p) where p = {k=1...len , index[k]<upperBound}*/
73
74
75 void Finley_Util_AddScatter(dim_t len,index_t* index,dim_t numData,double* in,double * out, index_t upperBound){
76 dim_t i,s;
77 for (s=0;s<len;s++) {
78 for(i=0;i<numData;i++) {
79 if( index[s]<upperBound ) {
80 out[INDEX2(i,index[s],numData)]+=in[INDEX2(i,s,numData)];
81 }
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 register double rtmp;
93 for (i=0;i<A1;i++) {
94 for (j=0;j<A2;j++) {
95 rtmp=0;
96 for (s=0;s<B2;s++) rtmp+=B[INDEX2(i,s,A1)]*C[INDEX2(s,j,B2)];
97 A[INDEX2(i,j,A1)]=rtmp;
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 register double rtmp;
109 for (q=0;q<len;q++) {
110 for (i=0;i<A1;i++) {
111 for (j=0;j<A2;j++) {
112 rtmp=0;
113 for (s=0;s<B2;s++) rtmp+=B[INDEX3(i,s,q,A1,B2)]*C[INDEX3(s,j,q,B2,A2)];
114 A[INDEX3(i,j,q, A1,A2)]=rtmp;
115 }
116 }
117 }
118 }
119 /* multiplies a set of matries with a single matrix: */
120
121 /* A(1:A1,1:A2,i)=B(1:A1,1:B2,i)*C(1:B2,1:A2) i=1,len */
122
123 void Finley_Util_SmallMatSetMult1(dim_t len,dim_t A1,dim_t A2, double* A, dim_t B2, double*B, double* C) {
124 dim_t q,i,j,s;
125 register double rtmp;
126 for (q=0;q<len;q++) {
127 for (i=0;i<A1;i++) {
128 for (j=0;j<A2;j++) {
129 rtmp=0;
130 for (s=0;s<B2;s++) rtmp+=B[INDEX3(i,s,q, A1,B2)]*C[INDEX2(s,j,B2)];
131 A[INDEX3(i,j,q,A1,A2)]=rtmp;
132 }
133 }
134 }
135 }
136 /* inverts the set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */
137 /* the determinante is returned. */
138
139 void Finley_Util_InvertSmallMat(dim_t len,dim_t dim,double* A,double *invA, double* det){
140 dim_t q;
141 register double D,A11,A12,A13,A21,A22,A23,A31,A32,A33;
142
143 switch(dim) {
144 case 1:
145 for (q=0;q<len;q++) {
146 D=A[q];
147 if (ABS(D) > 0 ){
148 det[q]=D;
149 D=1./D;
150 invA[q]=D;
151 } else {
152 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
153 return;
154 }
155 }
156 break;
157
158 case 2:
159 for (q=0;q<len;q++) {
160 A11=A[INDEX3(0,0,q,2,2)];
161 A12=A[INDEX3(0,1,q,2,2)];
162 A21=A[INDEX3(1,0,q,2,2)];
163 A22=A[INDEX3(1,1,q,2,2)];
164
165 D = A11*A22-A12*A21;
166 if (ABS(D) > 0 ){
167 det[q]=D;
168 D=1./D;
169 invA[INDEX3(0,0,q,2,2)]= A22*D;
170 invA[INDEX3(1,0,q,2,2)]=-A21*D;
171 invA[INDEX3(0,1,q,2,2)]=-A12*D;
172 invA[INDEX3(1,1,q,2,2)]= A11*D;
173 } else {
174 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
175 return;
176 }
177 }
178 break;
179
180 case 3:
181 for (q=0;q<len;q++) {
182 A11=A[INDEX3(0,0,q,3,3)];
183 A21=A[INDEX3(1,0,q,3,3)];
184 A31=A[INDEX3(2,0,q,3,3)];
185 A12=A[INDEX3(0,1,q,3,3)];
186 A22=A[INDEX3(1,1,q,3,3)];
187 A32=A[INDEX3(2,1,q,3,3)];
188 A13=A[INDEX3(0,2,q,3,3)];
189 A23=A[INDEX3(1,2,q,3,3)];
190 A33=A[INDEX3(2,2,q,3,3)];
191
192 D = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22);
193 if (ABS(D) > 0 ){
194 det[q] =D;
195 D=1./D;
196 invA[INDEX3(0,0,q,3,3)]=(A22*A33-A23*A32)*D;
197 invA[INDEX3(1,0,q,3,3)]=(A31*A23-A21*A33)*D;
198 invA[INDEX3(2,0,q,3,3)]=(A21*A32-A31*A22)*D;
199 invA[INDEX3(0,1,q,3,3)]=(A13*A32-A12*A33)*D;
200 invA[INDEX3(1,1,q,3,3)]=(A11*A33-A31*A13)*D;
201 invA[INDEX3(2,1,q,3,3)]=(A12*A31-A11*A32)*D;
202 invA[INDEX3(0,2,q,3,3)]=(A12*A23-A13*A22)*D;
203 invA[INDEX3(1,2,q,3,3)]=(A13*A21-A11*A23)*D;
204 invA[INDEX3(2,2,q,3,3)]=(A11*A22-A12*A21)*D;
205 } else {
206 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: Non-regular matrix");
207 return;
208 }
209 }
210 break;
211
212 }
213 return;
214 }
215
216 /* sets the derterminate of a set of dim x dim matrices A(:,:,1:len) with dim=1,2,3 */
217
218 void Finley_Util_DetOfSmallMat(dim_t len,dim_t dim,double* A, double* det){
219 dim_t q;
220 register double A11,A12,A13,A21,A22,A23,A31,A32,A33;
221
222 switch(dim) {
223 case 1:
224 for (q=0;q<len;q++) {
225 det[q]=A[q];
226 }
227 break;
228
229 case 2:
230 for (q=0;q<len;q++) {
231 A11=A[INDEX3(0,0,q,2,2)];
232 A12=A[INDEX3(0,1,q,2,2)];
233 A21=A[INDEX3(1,0,q,2,2)];
234 A22=A[INDEX3(1,1,q,2,2)];
235
236 det[q] = A11*A22-A12*A21;
237 }
238 break;
239
240 case 3:
241 for (q=0;q<len;q++) {
242 A11=A[INDEX3(0,0,q,3,3)];
243 A21=A[INDEX3(1,0,q,3,3)];
244 A31=A[INDEX3(2,0,q,3,3)];
245 A12=A[INDEX3(0,1,q,3,3)];
246 A22=A[INDEX3(1,1,q,3,3)];
247 A32=A[INDEX3(2,1,q,3,3)];
248 A13=A[INDEX3(0,2,q,3,3)];
249 A23=A[INDEX3(1,2,q,3,3)];
250 A33=A[INDEX3(2,2,q,3,3)];
251
252 det[q] = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22);
253 }
254 break;
255
256 }
257 return;
258 }
259 /* returns the normalized vector Normal[dim,len] orthogonal to A(:,0,q) and A(:,1,q) in the case of dim=3 */
260 /* or the vector A(:,0,q) in the case of dim=2 */
261
262 void Finley_NormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* Normal) {
263 dim_t q;
264 register double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33,length,invlength;
265
266 switch(dim) {
267 case 1:
268 for (q=0;q<len;q++) Normal[q] =1;
269 break;
270 case 2:
271 for (q=0;q<len;q++) {
272 A11=A[INDEX3(0,0,q,2,dim1)];
273 A21=A[INDEX3(1,0,q,2,dim1)];
274 length = sqrt(A11*A11+A21*A21);
275 if (! length>0) {
276 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero.");
277 return;
278 } else {
279 invlength=1./length;
280 Normal[INDEX2(0,q,2)]=A21*invlength;
281 Normal[INDEX2(1,q,2)]=-A11*invlength;
282 }
283 }
284 break;
285 case 3:
286 for (q=0;q<len;q++) {
287 A11=A[INDEX3(0,0,q,3,dim1)];
288 A21=A[INDEX3(1,0,q,3,dim1)];
289 A31=A[INDEX3(2,0,q,3,dim1)];
290 A12=A[INDEX3(0,1,q,3,dim1)];
291 A22=A[INDEX3(1,1,q,3,dim1)];
292 A32=A[INDEX3(2,1,q,3,dim1)];
293 CO_A13=A21*A32-A31*A22;
294 CO_A23=A31*A12-A11*A32;
295 CO_A33=A11*A22-A21*A12;
296 length=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33);
297 if (! length>0) {
298 Finley_setError(ZERO_DIVISION_ERROR,"__FILE__: area equals zero.");
299 return;
300 } else {
301 invlength=1./length;
302 Normal[INDEX2(0,q,3)]=CO_A13*invlength;
303 Normal[INDEX2(1,q,3)]=CO_A23*invlength;
304 Normal[INDEX2(2,q,3)]=CO_A33*invlength;
305 }
306
307 }
308 break;
309
310 }
311 return;
312 }
313
314 /* 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 */
315 /* or the vector A(:,0,q) in the case of dim=2 */
316
317 void Finley_LengthOfNormalVector(dim_t len, dim_t dim, dim_t dim1, double* A,double* length) {
318 dim_t q;
319 double A11,A12,CO_A13,A21,A22,CO_A23,A31,A32,CO_A33;
320
321 switch(dim) {
322 case 1:
323 for (q=0;q<len;q++) length[q] =1;
324 break;
325 case 2:
326 for (q=0;q<len;q++) {
327 A11=A[INDEX3(0,0,q,2,dim1)];
328 A21=A[INDEX3(1,0,q,2,dim1)];
329 length[q] = sqrt(A11*A11+A21*A21);
330 }
331 break;
332 case 3:
333 for (q=0;q<len;q++) {
334 A11=A[INDEX3(0,0,q,3,dim1)];
335 A21=A[INDEX3(1,0,q,3,dim1)];
336 A31=A[INDEX3(2,0,q,3,dim1)];
337 A12=A[INDEX3(0,1,q,3,dim1)];
338 A22=A[INDEX3(1,1,q,3,dim1)];
339 A32=A[INDEX3(2,1,q,3,dim1)];
340 CO_A13=A21*A32-A31*A22;
341 CO_A23=A31*A12-A11*A32;
342 CO_A33=A11*A22-A21*A12;
343 length[q]=sqrt(CO_A13*CO_A13+CO_A23*CO_A23+CO_A33*CO_A33);
344 }
345 break;
346
347 }
348 return;
349 }
350
351 /* inverts the map map of length len */
352 /* there is no range checking! */
353 /* at output Map[invMap[i]]=i for i=0:lenInvMap */
354
355 void Finley_Util_InvertMap(dim_t lenInvMap, index_t* invMap,dim_t lenMap, index_t* Map) {
356 dim_t i;
357 for (i=0;i<lenInvMap;i++) invMap[i]=0;
358 for (i=0;i<lenMap;i++) {
359 if (Map[i]>=0) invMap[Map[i]]=i;
360 }
361 }
362
363 /* orders a Finley_Util_ValueAndIndex array by value */
364 /* it is assumed that n is large */
365
366 int Finley_Util_ValueAndIndex_compar(const void *arg1 , const void *arg2 ) {
367 Finley_Util_ValueAndIndex *e1,*e2;
368 e1=(Finley_Util_ValueAndIndex*) arg1;
369 e2=(Finley_Util_ValueAndIndex*) arg2;
370 if (e1->value < e2->value) return -1;
371 if (e1->value > e2->value) return 1;
372 if (e1->index < e2->index) return -1;
373 if (e1->index > e2->index) return 1;
374 return 0;
375 }
376
377 void Finley_Util_sortValueAndIndex(dim_t n,Finley_Util_ValueAndIndex* array) {
378 /* OMP : needs parallelization !*/
379 qsort(array,n,sizeof(Finley_Util_ValueAndIndex),Finley_Util_ValueAndIndex_compar);
380 }
381
382
383 /**************************************************************/
384
385 /* calculates the minimum value from a dim X N integer array */
386
387 index_t Finley_Util_getMinInt(dim_t dim,dim_t N,index_t* values) {
388 dim_t i,j;
389 index_t out,out_local;
390 out=INDEX_T_MAX;
391 if (values!=NULL && dim*N>0 ) {
392 out=values[0];
393 #pragma omp parallel private(out_local)
394 {
395 out_local=out;
396 #pragma omp for private(i,j) schedule(static)
397 for (j=0;j<N;j++) {
398 for (i=0;i<dim;i++) out_local=MIN(out_local,values[INDEX2(i,j,dim)]);
399 }
400 #pragma omp critical
401 out=MIN(out_local,out);
402 }
403 }
404 return out;
405 }
406
407 /* calculates the maximum value from a dim X N integer array */
408
409 index_t Finley_Util_getMaxInt(dim_t dim,dim_t N,index_t* values) {
410 dim_t i,j;
411 index_t out,out_local;
412 out=-INDEX_T_MAX;
413 if (values!=NULL && dim*N>0 ) {
414 out=values[0];
415 #pragma omp parallel private(out_local)
416 {
417 out_local=out;
418 #pragma omp for private(i,j) schedule(static)
419 for (j=0;j<N;j++) {
420 for (i=0;i<dim;i++) out_local=MAX(out_local,values[INDEX2(i,j,dim)]);
421 }
422 #pragma omp critical
423 out=MAX(out_local,out);
424 }
425 }
426 return out;
427 }
428 /**************************************************************/
429
430 /* calculates the minimum value from a dim X N integer array */
431
432 index_t Finley_Util_getFlaggedMinInt(dim_t dim,dim_t N,index_t* values, index_t ignore) {
433 dim_t i,j;
434 index_t out,out_local;
435 out=INDEX_T_MAX;
436 if (values!=NULL && dim*N>0 ) {
437 out=values[0];
438 #pragma omp parallel private(out_local)
439 {
440 out_local=out;
441 #pragma omp for private(i,j) schedule(static)
442 for (j=0;j<N;j++) {
443 for (i=0;i<dim;i++) if (values[INDEX2(i,j,dim)]!=ignore) out_local=MIN(out_local,values[INDEX2(i,j,dim)]);
444 }
445 #pragma omp critical
446 out=MIN(out_local,out);
447 }
448 }
449 return out;
450 }
451
452 /* calculates the maximum value from a dim X N integer array */
453
454 index_t Finley_Util_getFlaggedMaxInt(dim_t dim,dim_t N,index_t* values, index_t ignore) {
455 dim_t i,j;
456 index_t out,out_local;
457 out=-INDEX_T_MAX;
458 if (values!=NULL && dim*N>0 ) {
459 out=values[0];
460 #pragma omp parallel private(out_local)
461 {
462 out_local=out;
463 #pragma omp for private(i,j) schedule(static)
464 for (j=0;j<N;j++) {
465 for (i=0;i<dim;i++) if (values[INDEX2(i,j,dim)]!=ignore) out_local=MAX(out_local,values[INDEX2(i,j,dim)]);
466 }
467 #pragma omp critical
468 out=MAX(out_local,out);
469 }
470 }
471 return out;
472 }
473
474 /* set the index of the positive entries in mask. The length of index is returned. */
475
476 dim_t Finley_Util_packMask(dim_t N,index_t* mask,index_t* index) {
477 dim_t out,k;
478 out=0;
479 /*OMP */
480 for (k=0;k<N;k++) {
481 if (mask[k]>=0) {
482 index[out]=k;
483 out++;
484 }
485 }
486 return out;
487 }
488
489 /* returns true if array contains value */
490 bool_t Finley_Util_isAny(dim_t N,index_t* array,index_t value) {
491 bool_t out=FALSE;
492 dim_t i;
493 #pragma omp parallel for private(i) schedule(static) reduction(||:out)
494 for (i=0;i<N;i++) out = out || (array[i]==value);
495 return out;
496 }
497 /* calculates the cummultative sum in array and returns the total sum */
498 index_t Finley_Util_cumsum(dim_t N,index_t* array) {
499 index_t out=0,tmp;
500 dim_t i;
501 #ifdef _OPENMP
502 index_t *partial_sums=NULL, sum;
503 partial_sums=TMPMEMALLOC(omp_get_max_threads(),index_t);
504 #pragma omp parallel private(sum,i,tmp)
505 {
506 sum=0;
507 #pragma omp for schedule(static)
508 for (i=0;i<N;++i) sum+=array[i];
509 partial_sums[omp_get_thread_num()]=sum;
510 #pragma omp barrier
511 #pragma omp master
512 {
513 out=0;
514 for (i=0;i<omp_get_max_threads();++i) {
515 tmp=out;
516 out+=partial_sums[i];
517 partial_sums[i]=tmp;
518 }
519 }
520 #pragma omp barrier
521 sum=partial_sums[omp_get_thread_num()];
522 #pragma omp for schedule(static)
523 for (i=0;i<N;++i) {
524 tmp=sum;
525 sum+=array[i];
526 array[i]=tmp;
527 }
528 }
529 TMPMEMFREE(partial_sums);
530 #else
531 for (i=0;i<N;++i) {
532 tmp=out;
533 out+=array[i];
534 array[i]=tmp;
535 }
536 #endif
537 return out;
538 }
539 void Finley_Util_setValuesInUse(const index_t *values, const dim_t numValues, dim_t *numValuesInUse, index_t **valuesInUse, Paso_MPIInfo* mpiinfo)
540 {
541 dim_t i;
542 index_t lastFoundValue=INDEX_T_MIN, minFoundValue, local_minFoundValue, *newValuesInUse=NULL;
543 register index_t itmp;
544 bool_t allFound=FALSE;
545 dim_t nv=0;
546
547 while (! allFound) {
548 /*
549 * find smallest value bigger than lastFoundValue
550 */
551 minFoundValue=INDEX_T_MAX;
552 #pragma omp parallel private(local_minFoundValue)
553 {
554 local_minFoundValue=minFoundValue;
555 #pragma omp for private(i,itmp) schedule(static)
556 for (i=0;i< numValues;i++) {
557 itmp=values[i];
558 if ((itmp>lastFoundValue) && (itmp<local_minFoundValue)) local_minFoundValue=itmp;
559 }
560 #pragma omp critical
561 {
562 if (local_minFoundValue<minFoundValue) minFoundValue=local_minFoundValue;
563 }
564
565 }
566 #ifdef PASO_MPI
567 local_minFoundValue=minFoundValue;
568 MPI_Allreduce(&local_minFoundValue,&minFoundValue, 1, MPI_INT, MPI_MIN, mpiinfo->comm );
569 #endif
570 /* if we found a new tag we need to add this too the valuesInUseList */
571
572 if (minFoundValue < INDEX_T_MAX) {
573 newValuesInUse=MEMALLOC(nv+1,index_t);
574 if (*valuesInUse!=NULL) {
575 memcpy(newValuesInUse,*valuesInUse,sizeof(index_t)*nv);
576 MEMFREE(*valuesInUse);
577 }
578 newValuesInUse[nv]=minFoundValue;
579 *valuesInUse=newValuesInUse;
580 newValuesInUse=NULL;
581 nv++;
582 lastFoundValue=minFoundValue;
583 } else {
584 allFound=TRUE;
585 }
586 }
587 *numValuesInUse=nv;
588 }
589
590
591 #ifdef PASO_MPI
592 void Finley_printDoubleArray( FILE *fid, dim_t n, double *array, char *name )
593 {
594 index_t i;
595
596 if( name )
597 fprintf( fid, "%s [ ", name );
598 else
599 fprintf( fid, "[ " );
600 for( i=0; i<(n<60 ? n : 60); i++ )
601 fprintf( fid, "%g ", array[i] );
602 if( n>=30 )
603 fprintf( fid, "... " );
604 fprintf( fid, "]\n" );
605 }
606 void Finley_printIntArray( FILE *fid, dim_t n, int *array, char *name )
607 {
608 index_t i;
609
610 if( name )
611 fprintf( fid, "%s [ ", name );
612 else
613 fprintf( fid, "[ " );
614 for( i=0; i<(n<60 ? n : 60); i++ )
615 fprintf( fid, "%d ", array[i] );
616 if( n>=30 )
617 fprintf( fid, "... " );
618 fprintf( fid, "]\n" );
619 }
620 void Finley_printMaskArray( FILE *fid, dim_t n, int *array, char *name )
621 {
622 index_t i;
623
624 if( name )
625 fprintf( fid, "%s [ ", name );
626 else
627 fprintf( fid, "[ " );
628 for( i=0; i<(n<60 ? n : 60); i++ )
629 if( array[i]!=-1 )
630 fprintf( fid, "%3d ", array[i] );
631 else
632 fprintf( fid, " * " );
633 if( n>=30 )
634 fprintf( fid, "... " );
635 fprintf( fid, "]\n" );
636 }
637 #endif

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26