/[escript]/branches/schroedinger_upto1946/finley/src/Assemble_getSize.c
ViewVC logotype

Contents of /branches/schroedinger_upto1946/finley/src/Assemble_getSize.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1947 - (show annotations)
Wed Oct 29 23:19:45 2008 UTC (10 years, 4 months ago) by jfenwick
File MIME type: text/plain
File size: 4985 byte(s)
This does not actually have the changes in it yet.

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 #pragma omp parallel private(local_X)
75 {
76 /* allocation of work arrays */
77 local_X=THREAD_MEMALLOC(NN*numDim,double);
78 if (! Finley_checkPtr(local_X) ) {
79 /* open the element loop */
80 #pragma omp for private(e,min_diff,diff,n0,n1,d,q,i,element_size_array) schedule(static)
81 for(e=0;e<elements->numElements;e++) {
82 /* gather local coordinates of nodes into local_X(numDim,NN): */
83 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
84 /* calculate minimal differences */
85 min_diff=-1;
86 for (n0=0;n0<NVertices;n0++) {
87 for (n1=n0+1;n1<NVertices;n1++) {
88 diff=0;
89 for (i=0;i<numDim;i++) {
90 d=local_X[INDEX2(i,n0,numDim)]-local_X[INDEX2(i,n1,numDim)];
91 diff+=d*d;
92 }
93 if (min_diff<0) {
94 min_diff=diff;
95 } else {
96 min_diff=MIN(min_diff,diff);
97 }
98 }
99 }
100 min_diff=sqrt(MAX(min_diff,0));
101 /* set all values to min_diff */
102 element_size_array=getSampleData(element_size,e);
103 for (q=0;q<numQuad;q++) element_size_array[q]=min_diff;
104 }
105 }
106 THREAD_MEMFREE(local_X);
107 }
108 }
109 return;
110 }
111 /*
112 * $Log$
113 * Revision 1.7 2005/09/15 03:44:21 jgs
114 * Merge of development branch dev-02 back to main trunk on 2005-09-15
115 *
116 * Revision 1.6.2.1 2005/09/07 06:26:17 gross
117 * the solver from finley are put into the standalone package paso now
118 *
119 * Revision 1.6 2005/07/08 04:07:47 jgs
120 * Merge of development branch back to main trunk on 2005-07-08
121 *
122 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
123 * some changes towards 64 integers in finley
124 *
125 * Revision 1.1.1.1.2.2 2005/03/02 23:35:05 gross
126 * reimplementation of the ILU in Finley. block size>1 still needs some testing
127 *
128 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
129 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
130 *
131 * Revision 1.1.1.1 2004/10/26 06:53:57 jgs
132 * initial import of project esys2
133 *
134 * Revision 1.3 2004/07/30 04:37:06 gross
135 * escript and finley are linking now and RecMeshTest.py has been passed
136 *
137 * Revision 1.2 2004/07/21 05:00:54 gross
138 * name changes in DataC
139 *
140 * Revision 1.1 2004/07/02 04:21:13 gross
141 * Finley C code has been included
142 *
143 *
144 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26