/[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 1072 - (show annotations)
Thu Mar 29 06:44:30 2007 UTC (12 years ago) by gross
File MIME type: text/plain
File size: 5325 byte(s)
PDE assemblage for reduced integration order + tests added.


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26