/[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 2595 - (show annotations)
Thu Aug 6 00:49:09 2009 UTC (9 years, 10 months ago) by gross
File MIME type: text/plain
File size: 3898 byte(s)
there was wrong test which triggered an excpetion in some cases where things were actually right
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 /* 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: Spatial and element dimension must match.");
64 } else */
65 if (! numSamplesEqual(element_size,numQuad,elements->numElements)) {
66 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: illegal number of samples of element size Data object");
67 } else if (! isDataPointShapeEqual(element_size,0,&(numDim))) {
68 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: illegal data point shape of element size Data object");
69 } else if (!isExpanded(element_size)) {
70 Finley_setError(TYPE_ERROR,"Finley_Assemble_getSize: expanded Data object is expected for element size.");
71 }
72 /* now we can start: */
73
74 if (Finley_noError()) {
75 requireWrite(element_size);
76 #pragma omp parallel private(local_X)
77 {
78 /* allocation of work arrays */
79 local_X=THREAD_MEMALLOC(NN*numDim,double);
80 if (! Finley_checkPtr(local_X) ) {
81 /* open the element loop */
82 #pragma omp for private(e,min_diff,diff,n0,n1,d,q,i,element_size_array) schedule(static)
83 for(e=0;e<elements->numElements;e++) {
84 /* gather local coordinates of nodes into local_X(numDim,NN): */
85 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
86 /* calculate minimal differences */
87 min_diff=-1;
88 for (n0=0;n0<NVertices;n0++) {
89 for (n1=n0+1;n1<NVertices;n1++) {
90 diff=0;
91 for (i=0;i<numDim;i++) {
92 d=local_X[INDEX2(i,n0,numDim)]-local_X[INDEX2(i,n1,numDim)];
93 diff+=d*d;
94 }
95 if (min_diff<0) {
96 min_diff=diff;
97 } else {
98 min_diff=MIN(min_diff,diff);
99 }
100 }
101 }
102 min_diff=sqrt(MAX(min_diff,0));
103 /* set all values to min_diff */
104 element_size_array=getSampleDataRW(element_size,e);
105 for (q=0;q<numQuad;q++) element_size_array[q]=min_diff;
106 }
107 }
108 THREAD_MEMFREE(local_X);
109 }
110 }
111 return;
112 }

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26