/[escript]/trunk/finley/src/Assemble_integrate.c
ViewVC logotype

Annotation of /trunk/finley/src/Assemble_integrate.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 155 - (hide annotations)
Wed Nov 9 02:02:19 2005 UTC (14 years ago) by jgs
Original Path: trunk/finley/src/finleyC/Assemble_integrate.c
File MIME type: text/plain
File size: 6565 byte(s)
move all directories from trunk/esys2 into trunk and remove esys2

1 jgs 150 /*
2     ******************************************************************************
3     * *
4     * COPYRIGHT ACcESS 2003,2004,2005 - All Rights Reserved *
5     * *
6     * This software is the property of ACcESS. No part of this code *
7     * may be copied in any form or by any means without the expressed written *
8     * consent of ACcESS. Copying, use or modification of this software *
9     * by any unauthorised person is illegal unless that person has a software *
10     * license agreement with ACcESS. *
11     * *
12     ******************************************************************************
13     */
14 jgs 82
15     /**************************************************************/
16    
17     /* assemblage routines: integrates data on quadrature points */
18    
19     /**************************************************************/
20    
21 jgs 150 /* Copyrights by ACcESS Australia 2003,2004,2005 */
22     /* Author: gross@access.edu.au */
23     /* version: $Id$ */
24 jgs 82
25     /**************************************************************/
26    
27     #include "Assemble.h"
28     #include "Util.h"
29     #ifdef _OPENMP
30     #include <omp.h>
31     #endif
32    
33     /**************************************************************/
34    
35     void Finley_Assemble_integrate(Finley_NodeFile* nodes, Finley_ElementFile* elements,escriptDataC* data,double* out) {
36 jgs 123 type_t datacase;
37     index_t node_offset;
38     dim_t q,e,i;
39 jgs 82 double *local_X, *dVdv, *Vol,*out_local,*data_array,rtmp;
40     if (nodes==NULL || elements==NULL) return;
41 jgs 123 dim_t NN=elements->ReferenceElement->Type->numNodes;
42     dim_t NS=elements->ReferenceElement->Type->numShapes;
43     type_t data_type=getFunctionSpaceType(data);
44     dim_t numComps=getDataPointSize(data);
45     dim_t numQuad=elements->ReferenceElement->numQuadNodes;
46     dim_t numDim=nodes->numDim;
47     dim_t numDim_local=elements->ReferenceElement->Type->numDim;
48 jgs 150 Finley_resetError();
49 jgs 82
50     /* set some parameter */
51    
52     if (data_type==FINLEY_ELEMENTS) {
53     datacase=0;
54     node_offset=0;
55     } else if (data_type==FINLEY_FACE_ELEMENTS) {
56     datacase=1;
57     node_offset=0;
58     } else if (data_type==FINLEY_CONTACT_ELEMENTS_1) {
59     datacase=1;
60     node_offset=0;
61     } else if (data_type==FINLEY_CONTACT_ELEMENTS_2) {
62     datacase=1;
63     node_offset=NN-NS;
64     } else {
65 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: integration of data is not possible.");
66 jgs 82 }
67    
68     /* check the shape of the data */
69    
70     if (! numSamplesEqual(data,numQuad,elements->numElements)) {
71 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of integrant kernel Data object");
72 jgs 82 }
73    
74     /* now we can start */
75    
76 jgs 150 if (Finley_noError()) {
77 jgs 82 for (q=0;q<numComps;q++) out[q]=0;
78     #pragma omp parallel private(local_X,dVdv,Vol,out_local,i)
79     {
80    
81     /* allocation of work arrays */
82    
83     local_X=dVdv=Vol=out_local=NULL;
84 jgs 102 out_local=THREAD_MEMALLOC(numComps,double);
85     dVdv=THREAD_MEMALLOC(numQuad*numDim_local*numDim,double);
86     Vol=THREAD_MEMALLOC(numQuad,double);
87     local_X=THREAD_MEMALLOC(NS*numDim,double);
88 jgs 82
89     if (! (Finley_checkPtr(Vol) || Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) || Finley_checkPtr(out_local))) {
90    
91     /* initialize local result */
92    
93     for (i=0;i<numComps;i++) out_local[i]=0;
94    
95     /* open the element loop */
96    
97     #pragma omp for private(e,q,data_array,rtmp) schedule(static)
98     for(e=0;e<elements->numElements;e++) {
99     /* gather local coordinates of nodes into local_X: */
100     Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
101     /* calculate dVdv(i,j,q)=V(i,k)*DSDv(k,j,q) */
102     Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
103     /* calculate volume /area elements */
104     switch (datacase) {
105     case 0:
106     Finley_Util_DetOfSmallMat(numQuad,numDim,dVdv,Vol);
107     break;
108     case 1:
109     Finley_LengthOfNormalVector(numQuad,numDim,numDim_local,dVdv,Vol);
110     break;
111     }
112     /* out=out+ integrals */
113     data_array=getSampleData(data,e);
114     if (isExpanded(data)) {
115     for (q=0;q<numQuad;q++) {
116     rtmp=ABS(Vol[q])*elements->ReferenceElement->QuadWeights[q];
117     for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*rtmp;
118     }
119     } else {
120     rtmp=0.;
121     for (q=0;q<numQuad;q++) rtmp+=ABS(Vol[q])*elements->ReferenceElement->QuadWeights[q];
122     for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
123     }
124     }
125    
126     /* add local results to global result */
127     #pragma omp critical
128     for (i=0;i<numComps;i++) out[i]+=out_local[i];
129    
130     }
131     THREAD_MEMFREE(local_X);
132     THREAD_MEMFREE(dVdv);
133     THREAD_MEMFREE(Vol);
134     THREAD_MEMFREE(out_local);
135     }
136     }
137     }
138    
139     /*
140     * $Log$
141 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
142     * Merge of development branch dev-02 back to main trunk on 2005-09-15
143     *
144     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
145     * the solver from finley are put into the standalone package paso now
146     *
147 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
148     * Merge of development branch back to main trunk on 2005-07-08
149     *
150 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
151 jgs 97 * *** empty log message ***
152 jgs 123 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
153     * some changes towards 64 integers in finley
154 jgs 82 *
155 jgs 123 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
156     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
157 jgs 97 *
158 jgs 82 *
159 jgs 123 *
160 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26