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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 123 - (hide annotations)
Fri Jul 8 04:08:13 2005 UTC (14 years, 3 months ago) by jgs
Original Path: trunk/esys2/finley/src/finleyC/Assemble_getSize.c
File MIME type: text/plain
File size: 4486 byte(s)
Merge of development branch back to main trunk on 2005-07-08

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26