/[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 1388 - (hide annotations)
Fri Jan 11 07:45:58 2008 UTC (11 years, 10 months ago) by trankine
File MIME type: text/plain
File size: 4482 byte(s)
And get the *(&(*&(* name right
1 jgs 82
2 ksteube 1312 /* $Id$ */
3 jgs 82
4 ksteube 1312 /*******************************************************
5     *
6     * Copyright 2003-2007 by ACceSS MNRF
7     * Copyright 2007 by University of Queensland
8     *
9     * http://esscc.uq.edu.au
10     * Primary Business: Queensland, Australia
11     * Licensed under the Open Software License version 3.0
12     * http://www.opensource.org/licenses/osl-3.0.php
13     *
14     *******************************************************/
15 jgs 82
16     /**************************************************************/
17    
18 ksteube 1312 /* assemblage routines: integrates data on quadrature points */
19 jgs 82
20     /**************************************************************/
21    
22     #include "Assemble.h"
23     #include "Util.h"
24     #ifdef _OPENMP
25     #include <omp.h>
26     #endif
27    
28     /**************************************************************/
29    
30     void Finley_Assemble_integrate(Finley_NodeFile* nodes, Finley_ElementFile* elements,escriptDataC* data,double* out) {
31 jgs 123 type_t datacase;
32     index_t node_offset;
33 gross 777 bool_t reducedIntegrationOrder;
34 jgs 123 type_t data_type=getFunctionSpaceType(data);
35     dim_t numComps=getDataPointSize(data);
36 gross 1028 Finley_ElementFile_Jacobeans* jac=NULL;
37 ksteube 1312 Paso_MPI_rank my_mpi_rank;
38    
39 jgs 150 Finley_resetError();
40 gross 1028 if (nodes==NULL || elements==NULL) return;
41 ksteube 1312 my_mpi_rank = nodes->MPIInfo->rank;
42 jgs 82 /* set some parameter */
43 gross 1064 jac=Finley_ElementFile_borrowJacobeans(elements,nodes,FALSE,Finley_Assemble_reducedIntegrationOrder(data));
44 gross 777 if (Finley_noError()) {
45 jgs 82
46 gross 777 /* check the shape of the data */
47     if (! numSamplesEqual(data,jac->ReferenceElement->numQuadNodes,elements->numElements)) {
48     Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: illegal number of samples of integrant kernel Data object");
49     }
50     /* now we can start */
51 jgs 82
52 gross 777 if (Finley_noError()) {
53 gross 853 dim_t q,e,i;
54     double *out_local=NULL, rtmp,*data_array=NULL;
55 gross 777 for (q=0;q<numComps;q++) out[q]=0;
56 gross 853 #pragma omp parallel private(q,i,rtmp,data_array,out_local)
57 gross 777 {
58 gross 853 out_local=THREAD_MEMALLOC(numComps,double);
59     if (! Finley_checkPtr(out_local) ) {
60     /* initialize local result */
61 jgs 82
62 gross 853 for (i=0;i<numComps;i++) out_local[i]=0;
63 jgs 82
64 gross 853 /* open the element loop */
65    
66     if (isExpanded(data)) {
67     #pragma omp for private(e) schedule(static)
68     for(e=0;e<elements->numElements;e++) {
69 ksteube 1312 if (elements->Owner[e] == my_mpi_rank) {
70 gross 853 data_array=getSampleData(data,e);
71     for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) {
72     for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
73     }
74 ksteube 1312 }
75 gross 853 }
76     } else {
77     #pragma omp for private(e) schedule(static)
78     for(e=0;e<elements->numElements;e++) {
79 ksteube 1312 if (elements->Owner[e] == my_mpi_rank) {
80 gross 853 data_array=getSampleData(data,e);
81     rtmp=0.;
82     for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) rtmp+=jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
83     for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
84 ksteube 1312 }
85 gross 853 }
86 gross 777 }
87 gross 853 /* add local results to global result */
88     #pragma omp critical
89     for (i=0;i<numComps;i++) out[i]+=out_local[i];
90 gross 777 }
91 gross 853 THREAD_MEMFREE(out_local);
92 gross 777 }
93 jgs 82 }
94 gross 777 }
95 jgs 82 }
96    
97     /*
98     * $Log$
99 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
100     * Merge of development branch dev-02 back to main trunk on 2005-09-15
101     *
102     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
103     * the solver from finley are put into the standalone package paso now
104     *
105 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
106     * Merge of development branch back to main trunk on 2005-07-08
107     *
108 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
109 jgs 97 * *** empty log message ***
110 jgs 123 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
111     * some changes towards 64 integers in finley
112 jgs 82 *
113 jgs 123 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
114     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
115 jgs 97 *
116 jgs 82 *
117 jgs 123 *
118 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26