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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 853 - (show annotations)
Wed Sep 20 05:56:36 2006 UTC (13 years ago) by gross
File MIME type: text/plain
File size: 5166 byte(s)
some performance problems wit openmp fixed.
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 /* assemblage routines: */
16
17 /* calculates the minimum distance between two vertices of elements and assigns the value to each */
18 /* quadrature point in element_size */
19
20
21 /**************************************************************/
22
23 /* Copyrights by ACcESS Australia, 2003,2004,2005 */
24 /* author: gross@access.edu.au */
25 /* version: $Id$ */
26
27 /**************************************************************/
28
29 #include "Assemble.h"
30 #include "Util.h"
31 #ifdef _OPENMP
32 #include <omp.h>
33 #endif
34
35 /**************************************************************/
36 void Finley_Assemble_getSize(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* element_size) {
37
38 double *local_X=NULL,*element_size_array;
39 dim_t e,n0,n1,q,i;
40 index_t node_offset;
41 double d,diff,min_diff;
42 Finley_resetError();
43
44 if (nodes==NULL || elements==NULL) return;
45 dim_t NVertices=elements->ReferenceElement->Type->numVertices;
46 dim_t NN=elements->ReferenceElement->Type->numNodes;
47 dim_t NS=elements->ReferenceElement->Type->numShapes;
48 dim_t numQuad=elements->ReferenceElement->numQuadNodes;
49 dim_t numDim=nodes->numDim;
50
51 /* set a few more parameters */
52
53 if (getFunctionSpaceType(element_size)==FINLEY_CONTACT_ELEMENTS_2) {
54 node_offset=NN-NS;
55 } else {
56 node_offset=0;
57 }
58
59 /* check the dimensions of element_size */
60
61 if (numDim!=elements->ReferenceElement->Type->numDim) {
62 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: Gradient: Spatial and element dimension must match.");
63 } else if (! numSamplesEqual(element_size,numQuad,elements->numElements)) {
64 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: illegal number of samples of element size Data object");
65 } else if (! isDataPointShapeEqual(element_size,0,&(numDim))) {
66 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: illegal data point shape of element size Data object");
67 } else if (!isExpanded(element_size)) {
68 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: expanded Data object is expected for element size.");
69 }
70 /* now we can start: */
71
72 if (Finley_noError()) {
73 #pragma omp parallel private(local_X)
74 {
75 /* allocation of work arrays */
76 local_X=THREAD_MEMALLOC(NN*numDim,double);
77 if (! Finley_checkPtr(local_X) ) {
78 /* open the element loop */
79 #pragma omp for private(e,min_diff,diff,n0,n1,d,q,i,element_size_array) schedule(static)
80 for(e=0;e<elements->numElements;e++) {
81 /* gather local coordinates of nodes into local_X(numDim,NN): */
82 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
83 /* calculate minimal differences */
84 min_diff=-1;
85 for (n0=0;n0<NVertices;n0++) {
86 for (n1=n0+1;n1<NVertices;n1++) {
87 diff=0;
88 for (i=0;i<numDim;i++) {
89 d=local_X[INDEX2(i,n0,numDim)]-local_X[INDEX2(i,n1,numDim)];
90 diff+=d*d;
91 }
92 if (min_diff<0) {
93 min_diff=diff;
94 } else {
95 min_diff=MIN(min_diff,diff);
96 }
97 }
98 }
99 min_diff=sqrt(MAX(min_diff,0));
100 /* set all values to min_diff */
101 element_size_array=getSampleData(element_size,e);
102 for (q=0;q<numQuad;q++) element_size_array[q]=min_diff;
103 }
104 }
105 THREAD_MEMFREE(local_X);
106 }
107 }
108 return;
109 }
110 /*
111 * $Log$
112 * Revision 1.7 2005/09/15 03:44:21 jgs
113 * Merge of development branch dev-02 back to main trunk on 2005-09-15
114 *
115 * Revision 1.6.2.1 2005/09/07 06:26:17 gross
116 * the solver from finley are put into the standalone package paso now
117 *
118 * Revision 1.6 2005/07/08 04:07:47 jgs
119 * Merge of development branch back to main trunk on 2005-07-08
120 *
121 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
122 * some changes towards 64 integers in finley
123 *
124 * Revision 1.1.1.1.2.2 2005/03/02 23:35:05 gross
125 * reimplementation of the ILU in Finley. block size>1 still needs some testing
126 *
127 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
128 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
129 *
130 * Revision 1.1.1.1 2004/10/26 06:53:57 jgs
131 * initial import of project esys2
132 *
133 * Revision 1.3 2004/07/30 04:37:06 gross
134 * escript and finley are linking now and RecMeshTest.py has been passed
135 *
136 * Revision 1.2 2004/07/21 05:00:54 gross
137 * name changes in DataC
138 *
139 * Revision 1.1 2004/07/02 04:21:13 gross
140 * Finley C code has been included
141 *
142 *
143 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26