source: codes/icosagcm/devel/Python/src/unstructured.pyx @ 674

Last change on this file since 674 was 674, checked in by dubos, 6 years ago

devel/unstructured : mesh reordering (Morton)

File size: 14.4 KB
Line 
1import time
2import math
3import dynamico.wrap as wrap
4from dynamico.libs import libicosa
5
6from ctypes import c_void_p, c_int, c_double, c_bool
7
8#from libc.time cimport time as ctime, time_t
9from numpy cimport ndarray
10cimport numpy as np
11
12#------------- direct Cython interface to DYNAMICO routines -------------#
13
14cdef enum: max_nb_stage=5
15cdef extern :
16    cdef double tauj[max_nb_stage]
17    cdef double cslj[max_nb_stage][max_nb_stage]
18    cdef double cflj[max_nb_stage][max_nb_stage]
19    cdef int nb_stage[1]
20
21cdef extern from "functions.h":
22     cdef void dynamico_ARK_step(int nstep,
23                                 double *mass_col, double *rhodz, double *theta_rhodz, 
24                                 double *u, double *geopot, double *w,
25                                 double *theta, double *ps, double *pk, double *hflux, double *qv,
26                                 double *dmass_col, double *drhodz, double *dtheta_rhodz,
27                                 double *du_fast, double *du_slow,
28                                 double *dPhi_fast, double *dPhi_slow, 
29                                 double *dW_fast, double *dW_slow)
30     cdef void dynamico_init_params()
31     cpdef void dynamico_setup_xios()
32     cpdef void dynamico_xios_set_timestep(double)
33     cpdef void dynamico_xios_update_calendar(int)
34
35#------------- import and wrap DYNAMICO routines -------------#
36
37ker=wrap.Struct() # store imported fun X as funs.X
38
39check_args = False # use True instead of False for debugging, probably with some overhead
40
41try:   
42    kernels = wrap.SharedLib(vars(ker), libicosa, check_args=check_args) 
43    setvar, setvars, getvar, getvars = kernels.setvar, kernels.setvars, kernels.getvar, kernels.getvars
44except OSError:
45    print """
46Unable to load shared library 'libicosa.so' !
47    """
48    raise
49
50# providing a full prototype enables type-checking when calling
51# if a number n is present in the prototype, the previous type is repeated n times
52kernels.import_funs([
53                     ['dynamico_setup_xios',None],
54                     ['dynamico_xios_set_timestep',c_double],
55                     ['dynamico_xios_update_calendar',c_int],
56                     ['dynamico_init_mesh',c_void_p,13],
57                     ['dynamico_init_metric', c_void_p,6],
58                     ['dynamico_init_hybrid', c_void_p,3],
59                     ['dynamico_caldyn_unstructured', c_double, c_void_p,20],
60                     ['dynamico_partition_graph', c_int,2, c_void_p,3, c_int, c_void_p],
61                     ['dynamico_morton_encode', c_int,c_void_p,4]
62                     ])
63
64# set/get global variables
65eta_mass,eta_lag=(1,2)
66thermo_theta,thermo_entropy,thermo_moist,thermo_boussinesq=(1,2,3,4)
67
68kernels.addvars(
69    c_bool,'hydrostatic','debug_hevi_solver',
70    c_int,'llm','nqdyn','primal_num','max_primal_deg',
71    'dual_num','max_dual_deg','edge_num','max_trisk_deg',
72    'caldyn_thermo','caldyn_eta','nb_threads',
73    c_double,'elapsed','g', 'ptop', 'cpp', 'cppv',
74    'Rd', 'Rv', 'preff', 'Treff', 'pbot', 'rho_bot', 'Phi_bot')
75
76elapsed=0.
77
78#------------------------ Extension type performing a full ARK time step ----------------------
79
80ctypedef double *dbl_ptr
81cdef dbl_ptr ptr1(double[:] data) : return &data[0]
82cdef dbl_ptr ptr2(double[:,:] data) : return &data[0,0]
83cdef dbl_ptr ptr3(double[:,:,:] data) : return &data[0,0,0]
84cdef dbl_ptr ptr4(double[:,:,:,:] data) : return &data[0,0,0,0]
85cdef dbl_ptr ptr(data):
86    n=data.ndim
87    if n==1 : return ptr1(data)
88    if n==2 : return ptr2(data)
89    if n==3 : return ptr3(data)
90    if n==4 : return ptr4(data)
91    if n>4: raise IndexError
92       
93cdef alloc(dbl_ptr *p, allocator, n=1):
94    data=allocator(n)
95    p[0]=ptr(data)
96    return data
97
98cdef check_ptr(name, dbl_ptr p, ndarray data):
99    if p != ptr(data) : print name, 'p <> ptr(data) !!'
100       
101cdef class Caldyn_step:
102    # number of time steps to do at each invocation of advance()
103    cdef int nstep
104    # pointer to allocated arrays
105    cdef dbl_ptr p_mass, p_theta_rhodz, p_u, p_geopot, p_W                   # prognostic
106    cdef dbl_ptr p_mass_col, p_dmass_col, p_ps, p_theta, p_pk, p_hflux, p_qv # diagnostic
107    cdef dbl_ptr p_drhodz, p_dtheta_rhodz, p_du_fast, p_du_slow                # tendencies
108    cdef dbl_ptr p_dPhi_fast, p_dPhi_slow, p_dW_fast, p_dW_slow                # tendencies
109    # allocated arrays, must remain referenced or segfault
110    cdef readonly ndarray mass, theta_rhodz, u, geopot, W
111    cdef readonly ndarray mass_col, dmass_col, ps, theta, pk, hflux, qv
112    cdef readonly ndarray drhodz, dtheta_rhodz, du_fast, du_slow
113    cdef readonly ndarray dPhi_fast, dPhi_slow, dW_fast, dW_slow
114
115    def __init__(self,mesh,time_scheme, nstep):
116        self.nstep=nstep
117        #        self.mesh=mesh
118        fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass
119        fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z
120        # collect coefficients of time scheme
121        cdef double[:]   tauj_ = time_scheme.tauj
122        cdef double[:,:] cslj_ = time_scheme.csjl
123        cdef double[:,:] cflj_ = time_scheme.cfjl
124        ns = time_scheme.nstage
125        nb_stage[0]=ns
126
127        cdef int i,j
128        for i in range(ns):
129            tauj[i]=tauj_[i]
130            for j in range(ns):
131                cslj[i][j]=cslj_[i,j]
132                cflj[i][j]=cflj_[i,j]
133        # allocate arrays, store pointers to avoid overhead when calling dynamico
134        #    prognostic/diagnostic
135        self.ps                       = alloc(&self.p_ps, fps) 
136        self.mass_col, self.dmass_col = alloc(&self.p_mass_col, fps), alloc(&self.p_dmass_col, fps,ns), 
137        self.mass, self.theta_rhodz   = alloc(&self.p_mass, fmass), alloc(&self.p_theta_rhodz, fmass),
138        self.theta, self.pk           = alloc(&self.p_theta, fmass), alloc(&self.p_pk, fmass), 
139        self.geopot, self.W           = alloc(&self.p_geopot, fw), alloc(&self.p_W, fw),
140        self.hflux, self.u            = alloc(&self.p_hflux, fu), alloc(&self.p_u, fu) 
141        self.qv                       = alloc(&self.p_qv,fz)
142        #    tendencies
143        self.drhodz, self.dtheta_rhodz  = alloc(&self.p_drhodz,fmass,ns), alloc(&self.p_dtheta_rhodz,fmass,ns) 
144        self.du_fast, self.du_slow      = alloc(&self.p_du_fast,fu,ns), alloc(&self.p_du_slow,fu,ns)
145        self.dPhi_fast, self.dPhi_slow  = alloc(&self.p_dPhi_fast,fw,ns), alloc(&self.p_dPhi_slow,fw,ns)
146        self.dW_fast, self.dW_slow      = alloc(&self.p_dW_fast,fw,ns), alloc(&self.p_dW_slow,fw,ns)
147    def next(self):
148        #        global elapsed
149        # time1=time.time()
150        dynamico_ARK_step(self.nstep,
151                          self.p_mass_col, self.p_mass, self.p_theta_rhodz,
152                          self.p_u, self.p_geopot, self.p_W,
153                          self.p_theta, self.p_ps, self.p_pk, self.p_hflux, self.p_qv,
154                          self.p_dmass_col, self.p_drhodz, self.p_dtheta_rhodz,
155                          self.p_du_fast, self.p_du_slow,
156                          self.p_dPhi_fast, self.p_dPhi_slow,
157                          self.p_dW_fast, self.p_dW_slow)
158        #time2=time.time()
159        #if time2>time1: elapsed=elapsed+time2-time1
160
161def caldyn_step_TRSW(mesh,time_scheme,nstep):
162    setvars(('hydrostatic','caldyn_thermo','caldyn_eta'),
163            (True,thermo_boussinesq,eta_lag))
164    dynamico_init_params()
165    return Caldyn_step(mesh,time_scheme, nstep)
166def caldyn_step_HPE(mesh,time_scheme,nstep, caldyn_thermo,caldyn_eta, thermo,BC,g):
167    setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
168             'g','ptop','Rd','cpp','preff','Treff'),
169             (True,caldyn_thermo,caldyn_eta,
170              g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0))
171    dynamico_init_params()
172    return Caldyn_step(mesh,time_scheme, nstep)
173def caldyn_step_NH(mesh,time_scheme,nstep, caldyn_thermo, caldyn_eta, thermo,BC,g):
174    setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
175             'g','ptop','Rd','cpp','preff','Treff','pbot','rho_bot'),
176             (False,caldyn_thermo,caldyn_eta,
177              g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0,
178              BC.pbot.max(), BC.rho_bot.max()))
179    dynamico_init_params()
180    return Caldyn_step(mesh,time_scheme, nstep)
181   
182#----------------------------- Base class for dynamics ------------------------
183
184class Caldyn:
185    def __init__(self,mesh):
186        self.mesh=mesh
187        fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass
188        fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z
189        self.ps, self.ms, self.dms       = fps(), fps(), fps()
190        self.s, self.hs, self.dhs        = ftheta(), ftheta(), ftheta()
191        self.pk, self.berni, self.geopot, self.hflux = fmass(),fmass(),fw(),fu()
192        self.qu, self.qv                 = fu(),fz()
193        self.fmass, self.ftheta, self.fu, self.fw = fmass, ftheta, fu, fw
194    def bwd_fast_slow(self, flow, tau):
195        global elapsed
196        time1=time.time()
197        flow,fast,slow = self._bwd_fast_slow_(flow,tau)
198        time2=time.time()
199        elapsed=elapsed+time2-time1
200        return flow,fast,slow
201
202# when calling caldyn_unstructured, arrays for tendencies must be re-created each time
203# to avoid overwriting in the same memory space when time scheme is multi-stage
204
205#-------------------------- Shallow-water dynamics ---------------------
206
207class Caldyn_RSW(Caldyn):
208    def __init__(self,mesh):
209        Caldyn.__init__(self,mesh)
210        setvars(('hydrostatic','caldyn_thermo','caldyn_eta'),
211                (True,thermo_boussinesq,eta_lag))
212        self.dhs = self.fmass()
213        dynamico_init_params()
214    def _bwd_fast_slow_(self, flow, tau):
215        h,u = flow
216        # h*s = h => uniform buoyancy s=1 => shallow-water
217        dh, du_slow, du_fast, hs, buf = self.fmass(), self.fu(), self.fu(), h.copy(), self.geopot
218        ker.dynamico_caldyn_unstructured(tau, self.ms, h, hs, u, self.geopot, buf,
219                  self.s, self.ps, self.pk, self.hflux, self.qv,
220                  self.dms, dh, self.dhs, du_fast, du_slow,
221                  buf, buf, buf, buf)
222        return (h,u), (0.,du_fast), (dh,du_slow)
223
224#----------------------------------- HPE ------------------------------------
225
226class Caldyn_HPE(Caldyn):
227    def __init__(self,caldyn_thermo,caldyn_eta, mesh,thermo,BC,g):
228        Caldyn.__init__(self,mesh)
229        setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
230                 'g','ptop','Rd','cpp','preff','Treff'),
231                (True,caldyn_thermo,caldyn_eta,
232                 g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0))
233        dynamico_init_params()
234    def _bwd_fast_slow_(self, flow, tau):
235        dm, dS, du_slow, du_fast, buf = self.fmass(), self.ftheta(), self.fu(), self.fu(), self.geopot
236        m,S,u = flow
237        ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, self.geopot, buf,
238                  self.s, self.ps, self.pk, self.hflux, self.qv,
239                  self.dms, dm, dS, du_fast, du_slow,
240                  buf, buf, buf, buf)
241        return (m,S,u), (0.,0.,du_fast), (dm,dS,du_slow)
242
243#----------------------------------- NH ------------------------------------
244
245class Caldyn_NH(Caldyn):
246    def __init__(self,caldyn_thermo,caldyn_eta, mesh,thermo,BC,g):
247        Caldyn.__init__(self,mesh)
248        setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
249                 'g','ptop','Rd','cpp','preff','Treff',
250                 'pbot','rho_bot'),
251                (False,caldyn_thermo,caldyn_eta,
252                 g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0,
253                 BC.pbot.max(), BC.rho_bot.max()))
254        dynamico_init_params()
255    def bwd_fast_slow(self, flow, tau):
256        ftheta, fmass, fu, fw = self.ftheta, self.fmass, self.fu, self.fw
257        dm, dS, du_slow, du_fast = fmass(), ftheta(), fu(), fu()
258        dPhi_slow, dPhi_fast, dW_slow, dW_fast = fw(), fw(), fw(), fw()
259        m,S,u,Phi,W = flow
260        ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, Phi, W,
261                  self.s, self.ps, self.pk, self.hflux, self.qv,
262                  self.dms, dm, dS, du_fast, du_slow,
263                  dPhi_fast, dPhi_slow, dW_fast, dW_slow)
264        return ((m,S,u,Phi,W), (0.,0.,du_fast,dPhi_fast,dW_fast), 
265                (dm,dS,du_slow,dPhi_slow,dW_slow))
266
267#------------------------ Copy mesh info to Fortran side -------------------
268
269def init_mesh(llm, nqdyn, edge_num, primal_num, dual_num,
270              max_trisk_deg, max_primal_deg, max_dual_deg,
271              primal_nb, primal_edge, primal_ne,
272              dual_nb,dual_edge,dual_ne,dual_vertex, 
273              left,right,down,up,trisk_deg,trisk,
274              Ai, Av, fv, le_de, Riv2, wee):
275    setvars( ('llm','nqdyn','edge_num','primal_num','dual_num',
276              'max_trisk_deg','max_primal_deg','max_dual_deg'),
277             (llm, nqdyn, edge_num, primal_num, dual_num, 
278              max_trisk_deg, max_primal_deg, max_dual_deg) )       
279    print('init_mesh ...')
280    ker.dynamico_init_mesh(primal_nb,primal_edge,primal_ne,
281              dual_nb,dual_edge,dual_ne,dual_vertex,
282              left,right,down,up,trisk_deg,trisk)
283    print ('...done')
284    print('init_metric ...')
285    ker.dynamico_init_metric(Ai,Av,fv,le_de,Riv2,wee)
286    print ('...done')
287
288#------------------------ Mesh partitioning ------------------------
289
290# Helper functions and interface to ParMETIS
291# list_stencil converts an adjacency graph from array format index[num_cells, MAX_EDGES] to compressed format
292# loc_stencil returns the start/end indices (vtxdist) expected by ParMETIS
293# i.e. index[start:end] with start=vtxdist[cell], end=vtxdist[cell+1] lists the edges of cell 'cell'
294
295import numpy as np
296
297def list_stencil(degree, stencil, cond=lambda x:True):
298    for i in range(degree.size):
299        for j in range(degree[i]):
300            s=stencil[i,j]
301            if cond(s): yield stencil[i,j]
302               
303def loc_stencil(degree, stencil):
304    loc=0
305    for i in range(degree.size):
306        yield loc
307        loc=loc+degree[i]
308    yield loc
309
310def partition_mesh(degree, stencil, nparts): 
311    # arguments : PArray1D and PArray2D describing mesh, number of desired partitions
312    dim_cell, degree, stencil = degree.dim, degree.data, stencil.data
313    comm, vtxdist, idx_start, idx_end = dim_cell.comm, dim_cell.vtxdist, dim_cell.start, dim_cell.end
314    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
315    adjncy_loc, xadj_loc = list_stencil(degree, stencil), loc_stencil(degree, stencil)
316    adjncy_loc, xadj_loc = [np.asarray(list(x), dtype=np.int32) for x in (adjncy_loc, xadj_loc)]
317    owner = np.zeros(idx_end-idx_start, dtype=np.int32);
318    ker.dynamico_partition_graph(mpi_rank, mpi_size, vtxdist, xadj_loc, adjncy_loc, nparts, owner)
319    return owner
Note: See TracBrowser for help on using the repository browser.