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

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

branch merged with trunk @1645. arch file (ep&mpi) added for ADA

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