source: XIOS/dev/dev_trunk_omp/extern/remap/src/libmapper.cpp @ 1619

Last change on this file since 1619 was 1619, checked in by yushan, 5 years ago

branch_openmp merged and NOT tested with trunk r1618

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