/[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 616 - (show annotations)
Wed Mar 22 02:46:56 2006 UTC (13 years, 7 months ago) by elspeth
File MIME type: text/plain
File size: 6243 byte(s)
Copyright added to more source files.

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 dim_t q,e,i;
37 double *local_X, *dVdv, *Vol,*out_local,*data_array,rtmp;
38 if (nodes==NULL || elements==NULL) return;
39 dim_t NN=elements->ReferenceElement->Type->numNodes;
40 dim_t NS=elements->ReferenceElement->Type->numShapes;
41 type_t data_type=getFunctionSpaceType(data);
42 dim_t numComps=getDataPointSize(data);
43 dim_t numQuad=elements->ReferenceElement->numQuadNodes;
44 dim_t numDim=nodes->numDim;
45 dim_t numDim_local=elements->ReferenceElement->Type->numDim;
46 Finley_resetError();
47
48 /* set some parameter */
49
50 if (data_type==FINLEY_ELEMENTS) {
51 datacase=0;
52 node_offset=0;
53 } else if (data_type==FINLEY_FACE_ELEMENTS) {
54 datacase=1;
55 node_offset=0;
56 } else if (data_type==FINLEY_CONTACT_ELEMENTS_1) {
57 datacase=1;
58 node_offset=0;
59 } else if (data_type==FINLEY_CONTACT_ELEMENTS_2) {
60 datacase=1;
61 node_offset=NN-NS;
62 } else {
63 Finley_setError(TYPE_ERROR,"__FILE__: integration of data is not possible.");
64 }
65
66 /* check the shape of the data */
67
68 if (! numSamplesEqual(data,numQuad,elements->numElements)) {
69 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of integrant kernel Data object");
70 }
71
72 /* now we can start */
73
74 if (Finley_noError()) {
75 for (q=0;q<numComps;q++) out[q]=0;
76 #pragma omp parallel private(local_X,dVdv,Vol,out_local,i)
77 {
78
79 /* allocation of work arrays */
80
81 local_X=dVdv=Vol=out_local=NULL;
82 out_local=THREAD_MEMALLOC(numComps,double);
83 dVdv=THREAD_MEMALLOC(numQuad*numDim_local*numDim,double);
84 Vol=THREAD_MEMALLOC(numQuad,double);
85 local_X=THREAD_MEMALLOC(NS*numDim,double);
86
87 if (! (Finley_checkPtr(Vol) || Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) || Finley_checkPtr(out_local))) {
88
89 /* initialize local result */
90
91 for (i=0;i<numComps;i++) out_local[i]=0;
92
93 /* open the element loop */
94
95 #pragma omp for private(e,q,data_array,rtmp) schedule(static)
96 for(e=0;e<elements->numElements;e++) {
97 /* gather local coordinates of nodes into local_X: */
98 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
99 /* calculate dVdv(i,j,q)=V(i,k)*DSDv(k,j,q) */
100 Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
101 /* calculate volume /area elements */
102 switch (datacase) {
103 case 0:
104 Finley_Util_DetOfSmallMat(numQuad,numDim,dVdv,Vol);
105 break;
106 case 1:
107 Finley_LengthOfNormalVector(numQuad,numDim,numDim_local,dVdv,Vol);
108 break;
109 }
110 /* out=out+ integrals */
111 data_array=getSampleData(data,e);
112 if (isExpanded(data)) {
113 for (q=0;q<numQuad;q++) {
114 rtmp=ABS(Vol[q])*elements->ReferenceElement->QuadWeights[q];
115 for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*rtmp;
116 }
117 } else {
118 rtmp=0.;
119 for (q=0;q<numQuad;q++) rtmp+=ABS(Vol[q])*elements->ReferenceElement->QuadWeights[q];
120 for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
121 }
122 }
123
124 /* add local results to global result */
125 #pragma omp critical
126 for (i=0;i<numComps;i++) out[i]+=out_local[i];
127
128 }
129 THREAD_MEMFREE(local_X);
130 THREAD_MEMFREE(dVdv);
131 THREAD_MEMFREE(Vol);
132 THREAD_MEMFREE(out_local);
133 }
134 }
135 }
136
137 /*
138 * $Log$
139 * Revision 1.6 2005/09/15 03:44:21 jgs
140 * Merge of development branch dev-02 back to main trunk on 2005-09-15
141 *
142 * Revision 1.5.2.1 2005/09/07 06:26:18 gross
143 * the solver from finley are put into the standalone package paso now
144 *
145 * Revision 1.5 2005/07/08 04:07:48 jgs
146 * Merge of development branch back to main trunk on 2005-07-08
147 *
148 * Revision 1.4 2004/12/15 07:08:32 jgs
149 * *** empty log message ***
150 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
151 * some changes towards 64 integers in finley
152 *
153 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
154 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
155 *
156 *
157 *
158 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26