/[escript]/branches/doubleplusgood/finley/src/Assemble_getSize.cpp
ViewVC logotype

Contents of /branches/doubleplusgood/finley/src/Assemble_getSize.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 4327 - (show annotations)
Wed Mar 20 05:09:11 2013 UTC (6 years, 6 months ago) by jfenwick
File size: 3935 byte(s)
some finley memory
1
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2013 by University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Open Software License version 3.0
9 * http://www.opensource.org/licenses/osl-3.0.php
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development since 2012 by School of Earth Sciences
13 *
14 *****************************************************************************/
15
16
17 /************************************************************************************/
18
19 /* assemblage routines: */
20
21 /* calculates the minimum distance between two vertices of elements and assigns the value to each */
22 /* quadrature point in element_size */
23
24
25 /************************************************************************************/
26
27 #include "Assemble.h"
28 #include "Util.h"
29 #ifdef _OPENMP
30 #include <omp.h>
31 #endif
32
33 /************************************************************************************/
34 void Finley_Assemble_getSize(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* element_size) {
35
36 Finley_ShapeFunction *shape=NULL;
37 Finley_ReferenceElement* refElement=NULL;
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,max_diff, f;
42 Finley_resetError();
43
44 if (nodes==NULL || elements==NULL) return;
45 refElement=Finley_ReferenceElementSet_borrowReferenceElement(elements->referenceElementSet,Finley_Assemble_reducedIntegrationOrder(element_size));
46 shape= refElement->Parametrization;
47 numDim=nodes->numDim;
48 numQuad=shape->numQuadNodes;
49 NN=elements->numNodes;
50 NS=refElement->Parametrization->Type->numShapes;
51 NVertices=refElement->Parametrization->Type->numVertices;
52
53
54 /* set a few more parameters */
55
56 if (getFunctionSpaceType(element_size)==FINLEY_CONTACT_ELEMENTS_2) {
57 node_offset=refElement->Type->offsets[1];
58 } else {
59 node_offset=refElement->Type->offsets[0];
60 }
61 f=pow(0.5, pow((double)(refElement->Type->numSubElements), 1./(double)(numDim))-1);
62
63 /* check the dimensions of element_size */
64
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=new double[NN*numDim];
80 if (! Finley_checkPtr(local_X) ) {
81 /* open the element loop */
82 #pragma omp for private(e,max_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 max_diff=0.;
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 max_diff=MAX(max_diff,diff);
96 }
97 }
98 max_diff=sqrt(max_diff)*f;
99 /* set all values to max_diff */
100 element_size_array=getSampleDataRW(element_size,e);
101 for (q=0;q<numQuad;q++) element_size_array[q]=max_diff;
102 }
103 }
104 delete [] local_X;
105 } /* end of parallel region */
106 }
107 return;
108 }

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26