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

Last change on this file since 802 was 802, checked in by dubos, 5 years ago

devel/unstructured : reduced, configurable log output

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