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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26