/[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 2271 - (show annotations)
Mon Feb 16 05:08:29 2009 UTC (10 years, 8 months ago) by jfenwick
File MIME type: text/plain
File size: 5016 byte(s)
Merging version 2269 to trunk

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26