source: XIOS/trunk/extern/remap/src/libmapper.cpp @ 694

Last change on this file since 694 was 694, checked in by rlacroix, 9 years ago

Fix compilation issues caused by the new "remap" library.

Use our MPI header so that C++ bindings are always disabled.

File size: 5.5 KB
Line 
1#include "mpi.hpp"
2#include <vector>
3#include <cassert>
4#include <cstring>
5#ifdef WRITE_TIMING
6#include <fstream>
7#endif
8
9#include "node.hpp"
10#include "grid.hpp"
11#include "circle.hpp" // cptRadius
12#include "elt.hpp"
13#include "meshutil.hpp" // cptArea
14#include "mapper.hpp"
15#include "cputime.hpp" // cputime
16
17using namespace sphereRemap ;
18
19/* mapper is a ponter to a global class instance whoes members are allocated in the first step (finding the sizes of the weight arrays)
20   and deallocated during the second step (computing the weights) */
21Mapper *mapper;
22
23
24/** xxx_bounds_yyy is of length n_vert_per_cell_xxx*n_cell_xxx
25   This function computes the weights and returns number of weights is returned through the last argument.
26   To get the actuall weights, use `remap_get_weights`.
27*/
28extern "C" void remap_get_num_weights(const double* src_bounds_lon, const double* src_bounds_lat, int n_vert_per_cell_src, int n_cell_src,
29                     const double* src_pole,
30                     const double* dst_bounds_lon, const double* dst_bounds_lat, int n_vert_per_cell_dst, int n_cell_dst,
31                     const double* dst_pole,
32                     int order, int* n_weights)
33{
34        assert(src_bounds_lon);
35        assert(src_bounds_lat);
36        assert(n_vert_per_cell_src >= 3);
37        assert(n_cell_src >= 4);
38        assert(dst_bounds_lon);
39        assert(dst_bounds_lat);
40        assert(n_vert_per_cell_dst >= 3);
41        assert(n_cell_dst >= 4);
42        assert(1 <= order && order <= 2);
43
44  mapper = new Mapper(MPI_COMM_WORLD);
45  mapper->setVerbosity(PROGRESS) ;
46  mapper->setSourceMesh(src_bounds_lon, src_bounds_lat, n_vert_per_cell_src, n_cell_src, src_pole ) ;
47  mapper->setTargetMesh(dst_bounds_lon, dst_bounds_lat, n_vert_per_cell_dst, n_cell_dst, dst_pole ) ;
48
49/*
50        srcGrid.pole = Coord(src_pole[0], src_pole[1], src_pole[2]);
51        tgtGrid.pole = Coord(dst_pole[0], dst_pole[1], dst_pole[2]);
52
53        int mpiRank, mpiSize;
54        MPI_Comm_rank(MPI_COMM_WORLD, &mpiRank);
55        MPI_Comm_size(MPI_COMM_WORLD, &mpiSize);
56
57        std::vector<Elt>  src_elt; src_elt.reserve(n_cell_src);
58        std::vector<Node> src_msh; src_msh.reserve(n_cell_src);
59        for (int i = 0; i < n_cell_src; i++)
60        {
61                int offs = i*n_vert_per_cell_src;
62                Elt elt(src_bounds_lon + offs, src_bounds_lat + offs, n_vert_per_cell_src);
63                elt.src_id.rank = mpiRank;
64                elt.src_id.ind = i;
65                src_elt.push_back(elt);
66                src_msh.push_back(Node(elt.x, cptRadius(elt), &src_elt.back()));
67                cptEltGeom(src_elt[i], Coord(src_pole[0], src_pole[1], src_pole[2]));
68        }
69        std::vector<Elt>  dst_elt; dst_elt.reserve(n_cell_dst);
70        std::vector<Node> dst_msh; dst_msh.reserve(n_cell_dst);
71        for (int i = 0; i < n_cell_dst; i++)
72        {
73                int offs = i*n_vert_per_cell_dst;
74                Elt elt(dst_bounds_lon + offs, dst_bounds_lat + offs, n_vert_per_cell_dst);
75                dst_elt.push_back(elt);
76                dst_msh.push_back(Node(elt.x, cptRadius(elt), &dst_elt.back()));
77                cptEltGeom(dst_elt[i], Coord(dst_pole[0], dst_pole[1], dst_pole[2]));
78        }
79        double tic = cputime();
80        mapper = new Mapper(MPI_COMM_WORLD);
81  mapper->setVerbosity(PROGRESS) ;
82        mapper->buildSSTree(src_msh, dst_msh);
83        double tac = cputime();
84        vector<double> timings = mapper->computeWeights(dst_elt, src_elt, order);
85*/
86
87  vector<double> timings = mapper->computeWeights(order);
88/*
89#ifdef WRITE_TIMING
90        int nGloWeights, gloSrcSize, gloDstSize;
91        int locSrcSize = src_elt.size();
92        int locDstSize = dst_elt.size();
93        MPI_Reduce(&(mapper->nWeights), &nGloWeights, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
94        MPI_Reduce(&locSrcSize, &gloSrcSize, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
95        MPI_Reduce(&locDstSize, &gloDstSize, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
96
97        if (mpiRank == 0)
98        {
99                ofstream timeout;
100                timeout.open("timings.txt", fstream::out | fstream::app);
101                timeout << mpiSize << " " << src_elt.size() << " " << dst_elt.size() << " " << nGloWeights
102                        << " " << tac-tic << " " << timings[0] << " " << timings[1] << " " << timings[2] << endl;
103                timeout.close();
104        }
105#endif
106*/
107
108        *n_weights = mapper->nWeights;
109
110}
111
112extern "C" void remap_get_barycentres_and_areas(const double* bounds_lon, const double* bounds_lat, int n_vert_per_cell, int n_cell,
113                     const double* pole,
114                     double* centre_lon, double* centre_lat, double* areas)
115{
116        for (int i = 0; i < n_cell; i++)
117        {
118                int offs = i*n_vert_per_cell;
119                Elt elt(bounds_lon + offs, bounds_lat + offs, n_vert_per_cell);
120                cptEltGeom(elt, Coord(pole[0], pole[1], pole[2]));
121                if (centre_lon && centre_lat) lonlat(elt.x, centre_lon[i], centre_lat[i]);
122                if (areas) areas[i] = elt.area;
123        }
124}
125
126/*
127extern "C" void remap_get_weights(double* weights, int* src_indices, int* src_ranks, int* dst_indices)
128{
129        memcpy(weights, mapper->remapMatrix, mapper->nWeights*sizeof(double));
130        memcpy(src_indices, mapper->srcAddress, mapper->nWeights*sizeof(int));
131        memcpy(src_ranks, mapper->srcRank, mapper->nWeights*sizeof(int));
132        memcpy(dst_indices, mapper->dstAddress, mapper->nWeights*sizeof(int));
133        delete mapper;
134}
135*/
136
137extern "C" void remap_get_weights(double* weights, int* src_indices, int* dst_indices)
138{
139        memcpy(weights, mapper->remapMatrix, mapper->nWeights*sizeof(double));
140        memcpy(src_indices, mapper->srcAddress, mapper->nWeights*sizeof(int));
141        memcpy(dst_indices, mapper->dstAddress, mapper->nWeights*sizeof(int));
142        delete mapper;
143}
144
145extern "C" void mpi_init()
146{
147        /*int argc = 0;
148        char **argv = NULL;
149        MPI_Init(&argc, &argv);*/
150        MPI_Init(NULL, NULL);
151}
152
153extern "C" int mpi_rank()
154{
155        int rank;
156        MPI_Comm_rank(MPI_COMM_WORLD, &rank);
157        return rank;
158}
159
160extern "C" int mpi_size()
161{
162        int size;
163        MPI_Comm_size(MPI_COMM_WORLD, &size);
164        return size;
165}
166
167extern "C" void mpi_finalize()
168{
169        MPI_Finalize();
170}
Note: See TracBrowser for help on using the repository browser.