/[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 1628 - (show annotations)
Fri Jul 11 13:12:46 2008 UTC (11 years, 2 months ago) by phornby
File MIME type: text/plain
File size: 4400 byte(s)

Merge in /branches/windows_from_1456_trunk_1620_merged_in branch.

You will find a preserved pre-merge trunk in tags under tags/trunk_at_1625.
That will be useful for diffing & checking on my stupidity.

Here is a list of the conflicts and their resolution at this
point in time.


=================================================================================
(LLWS == looks like white space).

finley/src/Assemble_addToSystemMatrix.c - resolve to branch - unused var. may be wrong.....
finley/src/CPPAdapter/SystemMatrixAdapter.cpp - resolve to branch - LLWS
finley/src/CPPAdapter/MeshAdapter.cpp - resolve to branch - LLWS
paso/src/PCG.c - resolve to branch - unused var fixes.
paso/src/SolverFCT.c - resolve to branch - LLWS
paso/src/FGMRES.c - resolve to branch - LLWS
paso/src/Common.h - resolve to trunk version. It's omp.h's include... not sure it's needed,
but for the sake of saftey.....
paso/src/Functions.c - resolve to branch version, indentation/tab removal and return error
on bad unimplemented Paso_FunctionCall.
paso/src/SolverFCT_solve.c - resolve to branch version, unused vars
paso/src/SparseMatrix_MatrixVector.c - resolve to branch version, unused vars.
escript/src/Utils.cpp - resloved to branch, needs WinSock2.h
escript/src/DataExpanded.cpp - resolved to branch version - LLWS
escript/src/DataFactory.cpp - resolve to branch version
=================================================================================

This currently passes tests on linux (debian), but is not checked on windows or Altix yet.

This checkin is to make a trunk I can check out for windows to do tests on it.

Known outstanding problem is in the operator=() method of exceptions
causing warning messages on the intel compilers.

May the God of doughnuts have mercy on my soul.


1
2 /* $Id$ */
3
4 /*******************************************************
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
16 /**************************************************************/
17
18 /* assemblage routines: integrates data on quadrature points */
19
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 type_t data_type=getFunctionSpaceType(data);
32 dim_t numComps=getDataPointSize(data);
33 Finley_ElementFile_Jacobeans* jac=NULL;
34 Paso_MPI_rank my_mpi_rank;
35
36 Finley_resetError();
37 if (nodes==NULL || elements==NULL) return;
38 my_mpi_rank = nodes->MPIInfo->rank;
39 /* set some parameter */
40 jac=Finley_ElementFile_borrowJacobeans(elements,nodes,FALSE,Finley_Assemble_reducedIntegrationOrder(data));
41 if (Finley_noError()) {
42
43 /* check the shape of the data */
44 if (! numSamplesEqual(data,jac->ReferenceElement->numQuadNodes,elements->numElements)) {
45 Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: illegal number of samples of integrant kernel Data object");
46 }
47 /* now we can start */
48
49 if (Finley_noError()) {
50 dim_t q,e,i;
51 double *out_local=NULL, rtmp,*data_array=NULL;
52 for (q=0;q<numComps;q++) out[q]=0;
53 #pragma omp parallel private(q,i,rtmp,data_array,out_local)
54 {
55 out_local=THREAD_MEMALLOC(numComps,double);
56 if (! Finley_checkPtr(out_local) ) {
57 /* initialize local result */
58
59 for (i=0;i<numComps;i++) out_local[i]=0;
60
61 /* open the element loop */
62
63 if (isExpanded(data)) {
64 #pragma omp for private(e) schedule(static)
65 for(e=0;e<elements->numElements;e++) {
66 if (elements->Owner[e] == my_mpi_rank) {
67 data_array=getSampleData(data,e);
68 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) {
69 for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
70 }
71 }
72 }
73 } else {
74 #pragma omp for private(e) schedule(static)
75 for(e=0;e<elements->numElements;e++) {
76 if (elements->Owner[e] == my_mpi_rank) {
77 data_array=getSampleData(data,e);
78 rtmp=0.;
79 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) rtmp+=jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
80 for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
81 }
82 }
83 }
84 /* add local results to global result */
85 #pragma omp critical
86 for (i=0;i<numComps;i++) out[i]+=out_local[i];
87 }
88 THREAD_MEMFREE(out_local);
89 }
90 }
91 }
92 }
93
94 /*
95 * $Log$
96 * Revision 1.6 2005/09/15 03:44:21 jgs
97 * Merge of development branch dev-02 back to main trunk on 2005-09-15
98 *
99 * Revision 1.5.2.1 2005/09/07 06:26:18 gross
100 * the solver from finley are put into the standalone package paso now
101 *
102 * Revision 1.5 2005/07/08 04:07:48 jgs
103 * Merge of development branch back to main trunk on 2005-07-08
104 *
105 * Revision 1.4 2004/12/15 07:08:32 jgs
106 * *** empty log message ***
107 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
108 * some changes towards 64 integers in finley
109 *
110 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
111 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
112 *
113 *
114 *
115 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26