/[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 853 - (hide annotations)
Wed Sep 20 05:56:36 2006 UTC (13 years, 4 months ago) by gross
File MIME type: text/plain
File size: 5285 byte(s)
some performance problems wit openmp fixed.
1 jgs 150 /*
2 elspeth 616 ************************************************************
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 jgs 150 */
12 jgs 82
13     /**************************************************************/
14    
15     /* assemblage routines: integrates data on quadrature points */
16    
17     /**************************************************************/
18    
19 jgs 150 /* Copyrights by ACcESS Australia 2003,2004,2005 */
20     /* Author: gross@access.edu.au */
21     /* version: $Id$ */
22 jgs 82
23     /**************************************************************/
24    
25     #include "Assemble.h"
26     #include "Util.h"
27     #ifdef _OPENMP
28     #include <omp.h>
29     #endif
30    
31     /**************************************************************/
32    
33     void Finley_Assemble_integrate(Finley_NodeFile* nodes, Finley_ElementFile* elements,escriptDataC* data,double* out) {
34 jgs 123 type_t datacase;
35     index_t node_offset;
36 gross 777 bool_t reducedIntegrationOrder;
37 jgs 82 if (nodes==NULL || elements==NULL) return;
38 jgs 123 type_t data_type=getFunctionSpaceType(data);
39     dim_t numComps=getDataPointSize(data);
40 jgs 150 Finley_resetError();
41 jgs 82
42     /* set some parameter */
43    
44     if (data_type==FINLEY_ELEMENTS) {
45 gross 777 reducedIntegrationOrder=FALSE;
46 jgs 82 } else if (data_type==FINLEY_FACE_ELEMENTS) {
47 gross 777 reducedIntegrationOrder=FALSE;
48 jgs 82 } else if (data_type==FINLEY_CONTACT_ELEMENTS_1) {
49 gross 777 reducedIntegrationOrder=FALSE;
50 jgs 82 } else if (data_type==FINLEY_CONTACT_ELEMENTS_2) {
51 gross 777 reducedIntegrationOrder=FALSE;
52 jgs 82 } else {
53 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: integration of data is not possible.");
54 jgs 82 }
55    
56 gross 777 /* get access to jacobean */
57     Finley_ElementFile_Jacobeans* jac=Finley_ElementFile_borrowJacobeans(elements,nodes,FALSE,reducedIntegrationOrder);
58 jgs 82
59 gross 777 if (Finley_noError()) {
60 jgs 82
61 gross 777 /* check the shape of the data */
62     if (! numSamplesEqual(data,jac->ReferenceElement->numQuadNodes,elements->numElements)) {
63     Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: illegal number of samples of integrant kernel Data object");
64     }
65 jgs 82
66    
67 gross 777 /* now we can start */
68 jgs 82
69 gross 777 if (Finley_noError()) {
70 gross 853 dim_t q,e,i;
71     double *out_local=NULL, rtmp,*data_array=NULL;
72 gross 777 for (q=0;q<numComps;q++) out[q]=0;
73 gross 853 #pragma omp parallel private(q,i,rtmp,data_array,out_local)
74 gross 777 {
75 gross 853 out_local=THREAD_MEMALLOC(numComps,double);
76     if (! Finley_checkPtr(out_local) ) {
77     /* initialize local result */
78 jgs 82
79 gross 853 for (i=0;i<numComps;i++) out_local[i]=0;
80 jgs 82
81 gross 853 /* open the element loop */
82    
83     if (isExpanded(data)) {
84     #pragma omp for private(e) schedule(static)
85     for(e=0;e<elements->numElements;e++) {
86     data_array=getSampleData(data,e);
87     for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) {
88     for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
89     }
90     }
91     } else {
92     #pragma omp for private(e) schedule(static)
93     for(e=0;e<elements->numElements;e++) {
94     data_array=getSampleData(data,e);
95     rtmp=0.;
96     for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) rtmp+=jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
97     for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
98     }
99 gross 777 }
100 gross 853 /* add local results to global result */
101     #pragma omp critical
102     for (i=0;i<numComps;i++) out[i]+=out_local[i];
103 gross 777 }
104 gross 853 THREAD_MEMFREE(out_local);
105 gross 777 }
106 jgs 82 }
107 gross 777 }
108 jgs 82 }
109    
110     /*
111     * $Log$
112 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
113     * Merge of development branch dev-02 back to main trunk on 2005-09-15
114     *
115     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
116     * the solver from finley are put into the standalone package paso now
117     *
118 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
119     * Merge of development branch back to main trunk on 2005-07-08
120     *
121 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
122 jgs 97 * *** empty log message ***
123 jgs 123 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
124     * some changes towards 64 integers in finley
125 jgs 82 *
126 jgs 123 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
127     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
128 jgs 97 *
129 jgs 82 *
130 jgs 123 *
131 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26