xref: /libCEED/python/ceed.py (revision 477729cfdbe3383ce841609affbc919aa37fc70d)
139b2de37Sjeremylt# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
239b2de37Sjeremylt# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
339b2de37Sjeremylt# reserved. See files LICENSE and NOTICE for details.
439b2de37Sjeremylt#
539b2de37Sjeremylt# This file is part of CEED, a collection of benchmarks, miniapps, software
639b2de37Sjeremylt# libraries and APIs for efficient high-order finite element and spectral
739b2de37Sjeremylt# element discretizations for exascale applications. For more information and
839b2de37Sjeremylt# source code availability see http://github.com/ceed.
939b2de37Sjeremylt#
1039b2de37Sjeremylt# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
1139b2de37Sjeremylt# a collaborative effort of two U.S. Department of Energy organizations (Office
1239b2de37Sjeremylt# of Science and the National Nuclear Security Administration) responsible for
1339b2de37Sjeremylt# the planning and preparation of a capable exascale ecosystem, including
1439b2de37Sjeremylt# software, applications, hardware, advanced system engineering and early
1539b2de37Sjeremylt# testbed platforms, in support of the nation's exascale computing imperative.
1639b2de37Sjeremylt
1739b2de37Sjeremyltfrom _ceed_cffi import ffi, lib
1839b2de37Sjeremyltimport sys
19*477729cfSJeremy L Thompsonimport os
2039b2de37Sjeremyltimport io
210a0da059Sjeremyltimport tempfile
2239b2de37Sjeremyltfrom abc import ABC
2339b2de37Sjeremyltfrom .ceed_vector import Vector
2439b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1
2569a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction
2639b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction
2739b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator
2839b2de37Sjeremyltfrom .ceed_constants import *
2939b2de37Sjeremylt
3039b2de37Sjeremylt# ------------------------------------------------------------------------------
317a7b0fa3SJed Brown
327a7b0fa3SJed Brown
3339b2de37Sjeremyltclass Ceed():
3439b2de37Sjeremylt    """Ceed: core components."""
3539b2de37Sjeremylt    # Attributes
3639b2de37Sjeremylt    _pointer = ffi.NULL
3739b2de37Sjeremylt
3839b2de37Sjeremylt    # Constructor
3939b2de37Sjeremylt    def __init__(self, resource="/cpu/self"):
4039b2de37Sjeremylt        # libCEED object
4139b2de37Sjeremylt        self._pointer = ffi.new("Ceed *")
4239b2de37Sjeremylt
4339b2de37Sjeremylt        # libCEED call
4439b2de37Sjeremylt        resourceAscii = ffi.new("char[]", resource.encode("ascii"))
45*477729cfSJeremy L Thompson        os.environ["CEED_ERROR_HANDLER"] = "return"
46*477729cfSJeremy L Thompson        err_code = lib.CeedInit(resourceAscii, self._pointer)
47*477729cfSJeremy L Thompson        if err_code:
48*477729cfSJeremy L Thompson            raise Exception("Error initializing backend resource: " + resource)
49*477729cfSJeremy L Thompson        lib.CeedSetErrorHandler(
50*477729cfSJeremy L Thompson            self._pointer[0], ffi.addressof(
51*477729cfSJeremy L Thompson                lib, "CeedErrorStore"))
5239b2de37Sjeremylt
5339b2de37Sjeremylt    # Representation
5439b2de37Sjeremylt    def __repr__(self):
5539b2de37Sjeremylt        return "<Ceed instance at " + hex(id(self)) + ">"
5639b2de37Sjeremylt
570a0da059Sjeremylt    # String conversion for print() to stdout
580a0da059Sjeremylt    def __str__(self):
590a0da059Sjeremylt        """View a Ceed via print()."""
600a0da059Sjeremylt
610a0da059Sjeremylt        # libCEED call
620a0da059Sjeremylt        with tempfile.NamedTemporaryFile() as key_file:
630a0da059Sjeremylt            with open(key_file.name, 'r+') as stream_file:
640a0da059Sjeremylt                stream = ffi.cast("FILE *", stream_file)
650a0da059Sjeremylt
66*477729cfSJeremy L Thompson                err_code = lib.CeedView(self._pointer[0], stream)
67*477729cfSJeremy L Thompson                self._check_error(err_code)
680a0da059Sjeremylt
690a0da059Sjeremylt                stream_file.seek(0)
700a0da059Sjeremylt                out_string = stream_file.read()
710a0da059Sjeremylt
720a0da059Sjeremylt        return out_string
730a0da059Sjeremylt
74*477729cfSJeremy L Thompson    # Error handler
75*477729cfSJeremy L Thompson    def _check_error(self, err_code):
76*477729cfSJeremy L Thompson        """Check return code and retrieve error message for non-zero code"""
77*477729cfSJeremy L Thompson        if (err_code):
78*477729cfSJeremy L Thompson            message = ffi.new("char **")
79*477729cfSJeremy L Thompson            lib.CeedGetErrorMessage(self._pointer[0], message)
80*477729cfSJeremy L Thompson            raise Exception(ffi.string(message[0]).decode("UTF-8"))
81*477729cfSJeremy L Thompson
8239b2de37Sjeremylt    # Get Resource
8339b2de37Sjeremylt    def get_resource(self):
8439b2de37Sjeremylt        """Get the full resource name for a Ceed context.
8539b2de37Sjeremylt
8639b2de37Sjeremylt           Returns:
8739b2de37Sjeremylt             resource: resource name"""
8839b2de37Sjeremylt
8939b2de37Sjeremylt        # libCEED call
9039b2de37Sjeremylt        resource = ffi.new("char **")
91*477729cfSJeremy L Thompson        err_code = lib.CeedGetResource(self._pointer[0], resource)
92*477729cfSJeremy L Thompson        self._check_error(err_code)
9339b2de37Sjeremylt
9439b2de37Sjeremylt        return ffi.string(resource[0]).decode("UTF-8")
9539b2de37Sjeremylt
9639b2de37Sjeremylt    # Preferred MemType
9739b2de37Sjeremylt    def get_preferred_memtype(self):
9839b2de37Sjeremylt        """Return Ceed preferred memory type.
9939b2de37Sjeremylt
10039b2de37Sjeremylt           Returns:
10139b2de37Sjeremylt             memtype: Ceed preferred memory type"""
10239b2de37Sjeremylt
10339b2de37Sjeremylt        # libCEED call
10439b2de37Sjeremylt        memtype = ffi.new("CeedMemType *", MEM_HOST)
105*477729cfSJeremy L Thompson        err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype)
106*477729cfSJeremy L Thompson        self._check_error(err_code)
10739b2de37Sjeremylt
10839b2de37Sjeremylt        return memtype[0]
10939b2de37Sjeremylt
11039b2de37Sjeremylt    # CeedVector
11139b2de37Sjeremylt    def Vector(self, size):
11239b2de37Sjeremylt        """Ceed Vector: storing and manipulating vectors.
11339b2de37Sjeremylt
11439b2de37Sjeremylt           Args:
11539b2de37Sjeremylt             size: length of vector
11639b2de37Sjeremylt
11739b2de37Sjeremylt           Returns:
11839b2de37Sjeremylt             vector: Ceed Vector"""
11939b2de37Sjeremylt
12039b2de37Sjeremylt        return Vector(self, size)
12139b2de37Sjeremylt
12239b2de37Sjeremylt    # CeedElemRestriction
123d979a051Sjeremylt    def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets,
124d979a051Sjeremylt                        memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
12539b2de37Sjeremylt        """Ceed ElemRestriction: restriction from local vectors to elements.
12639b2de37Sjeremylt
12739b2de37Sjeremylt           Args:
128d979a051Sjeremylt             nelem: number of elements described by the restriction
12939b2de37Sjeremylt             elemsize: size (number of nodes) per element
130d979a051Sjeremylt             ncomp: number of field components per interpolation node
131d979a051Sjeremylt                      (1 for scalar fields)
132d979a051Sjeremylt             compstride: Stride between components for the same L-vector "node".
133d979a051Sjeremylt                           Data for node i, component k can be found in the
134d979a051Sjeremylt                           L-vector at index [offsets[i] + k*compstride].
135d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
136d979a051Sjeremylt                       the elements and fields given by this restriction.
137d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
138d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
13939b2de37Sjeremylt                         for the unknowns corresponding to element i, where
140d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
141d979a051Sjeremylt                         [0, lsize - 1].
142d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
143d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
14439b2de37Sjeremylt
14539b2de37Sjeremylt           Returns:
14639b2de37Sjeremylt             elemrestriction: Ceed ElemRestiction"""
14739b2de37Sjeremylt
148d979a051Sjeremylt        return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
149d979a051Sjeremylt                               offsets, memtype=memtype, cmode=cmode)
15039b2de37Sjeremylt
151d979a051Sjeremylt    def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides):
152d979a051Sjeremylt        """Ceed Identity ElemRestriction: strided restriction from local vectors
153d979a051Sjeremylt             to elements.
15439b2de37Sjeremylt
15539b2de37Sjeremylt           Args:
156d979a051Sjeremylt             nelem: number of elements described by the restriction
15739b2de37Sjeremylt             elemsize: size (number of nodes) per element
158a8d32208Sjeremylt             ncomp: number of field components per interpolation node
159a8d32208Sjeremylt                      (1 for scalar fields)
160d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
161d979a051Sjeremylt                      the elements and fields given by this restriction.
16269a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
16369a53589Sjeremylt                         The data for node i, component j, element k in the
16469a53589Sjeremylt                         L-vector is given by
16569a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
16639b2de37Sjeremylt
16739b2de37Sjeremylt           Returns:
16869a53589Sjeremylt             elemrestriction: Ceed Strided ElemRestiction"""
16939b2de37Sjeremylt
1707a7b0fa3SJed Brown        return StridedElemRestriction(
171d979a051Sjeremylt            self, nelem, elemsize, ncomp, lsize, strides)
17239b2de37Sjeremylt
173d979a051Sjeremylt    def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride,
174d979a051Sjeremylt                               lsize, offsets, memtype=lib.CEED_MEM_HOST,
175d979a051Sjeremylt                               cmode=lib.CEED_COPY_VALUES):
176d979a051Sjeremylt        """Ceed Blocked ElemRestriction: blocked restriction from local vectors
177d979a051Sjeremylt             to elements.
17839b2de37Sjeremylt
17939b2de37Sjeremylt           Args:
180d979a051Sjeremylt             nelem: number of elements described by the restriction
18139b2de37Sjeremylt             elemsize: size (number of nodes) per element
18239b2de37Sjeremylt             blksize: number of elements in a block
183d979a051Sjeremylt             ncomp: number of field components per interpolation node
184d979a051Sjeremylt                       (1 for scalar fields)
185d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
186d979a051Sjeremylt                      the elements and fields given by this restriction.
187d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
188d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
18939b2de37Sjeremylt                         for the unknowns corresponding to element i, where
190d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
191d979a051Sjeremylt                         [0, lsize - 1]. The backend will permute and pad this
19239b2de37Sjeremylt                         array to the desired ordering for the blocksize, which is
19339b2de37Sjeremylt                         typically given by the backend. The default reordering is
19439b2de37Sjeremylt                         to interlace elements.
195d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
196d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
19739b2de37Sjeremylt
19839b2de37Sjeremylt           Returns:
19939b2de37Sjeremylt             elemrestriction: Ceed Blocked ElemRestiction"""
20039b2de37Sjeremylt
201d979a051Sjeremylt        return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp,
202d979a051Sjeremylt                                      compstride, lsize, offsets,
203d979a051Sjeremylt                                      memtype=memtype, cmode=cmode)
20439b2de37Sjeremylt
205d979a051Sjeremylt    def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp,
206d979a051Sjeremylt                                      lsize, strides):
207d979a051Sjeremylt        """Ceed Blocked Strided ElemRestriction: blocked and strided restriction
208d979a051Sjeremylt             from local vectors to elements.
20969a53589Sjeremylt
21069a53589Sjeremylt           Args:
211d979a051Sjeremylt             nelem: number of elements described in the restriction
21269a53589Sjeremylt             elemsize: size (number of nodes) per element
21369a53589Sjeremylt             blksize: number of elements in a block
21469a53589Sjeremylt             ncomp: number of field components per interpolation node
21569a53589Sjeremylt                      (1 for scalar fields)
216d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
217d979a051Sjeremylt                      the elements and fields given by this restriction.
21869a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
21969a53589Sjeremylt                         The data for node i, component j, element k in the
22069a53589Sjeremylt                         L-vector is given by
22169a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
22269a53589Sjeremylt
22369a53589Sjeremylt           Returns:
22469a53589Sjeremylt             elemrestriction: Ceed Strided ElemRestiction"""
22569a53589Sjeremylt
22669a53589Sjeremylt        return BlockedStridedElemRestriction(self, nelem, elemsize, blksize,
227d979a051Sjeremylt                                             ncomp, lsize, strides)
22869a53589Sjeremylt
22939b2de37Sjeremylt    # CeedBasis
23039b2de37Sjeremylt    def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
23139b2de37Sjeremylt                      qref1d, qweight1d):
23239b2de37Sjeremylt        """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
23339b2de37Sjeremylt             H^1 discretizations.
23439b2de37Sjeremylt
23539b2de37Sjeremylt           Args:
23639b2de37Sjeremylt             dim: topological dimension
23739b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
23839b2de37Sjeremylt             P1d: number of nodes in one dimension
23939b2de37Sjeremylt             Q1d: number of quadrature points in one dimension
240d979a051Sjeremylt             *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix
241d979a051Sjeremylt                          expressing the values of nodal basis functions at
242d979a051Sjeremylt                          quadrature points
243d979a051Sjeremylt             *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix
244d979a051Sjeremylt                        expressing the derivatives of nodal basis functions at
245d979a051Sjeremylt                        quadrature points
246d979a051Sjeremylt             *qref1d: Numpy array of length Q1d holding the locations of
247d979a051Sjeremylt                        quadrature points on the 1D reference element [-1, 1]
248d979a051Sjeremylt             *qweight1d: Numpy array of length Q1d holding the quadrature
249d979a051Sjeremylt                           weights on the reference element
25039b2de37Sjeremylt
25139b2de37Sjeremylt           Returns:
25239b2de37Sjeremylt             basis: Ceed Basis"""
25339b2de37Sjeremylt
25439b2de37Sjeremylt        return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
25539b2de37Sjeremylt                             qref1d, qweight1d)
25639b2de37Sjeremylt
25739b2de37Sjeremylt    def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode):
258d979a051Sjeremylt        """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange
259d979a051Sjeremylt             basis objects for H^1 discretizations.
26039b2de37Sjeremylt
26139b2de37Sjeremylt           Args:
26239b2de37Sjeremylt             dim: topological dimension
26339b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
26439b2de37Sjeremylt             P: number of Gauss-Lobatto nodes in one dimension.  The
26539b2de37Sjeremylt                  polynomial degree of the resulting Q_k element is k=P-1.
26639b2de37Sjeremylt             Q: number of quadrature points in one dimension
26739b2de37Sjeremylt             qmode: distribution of the Q quadrature points (affects order of
26839b2de37Sjeremylt                      accuracy for the quadrature)
26939b2de37Sjeremylt
27039b2de37Sjeremylt           Returns:
27139b2de37Sjeremylt             basis: Ceed Basis"""
27239b2de37Sjeremylt
27339b2de37Sjeremylt        return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode)
27439b2de37Sjeremylt
27539b2de37Sjeremylt    def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight):
276d979a051Sjeremylt        """Ceed H1 Basis: finite element non tensor-product basis for H^1
277d979a051Sjeremylt             discretizations.
27839b2de37Sjeremylt
27939b2de37Sjeremylt           Args:
28039b2de37Sjeremylt             topo: topology of the element, e.g. hypercube, simplex, etc
28139b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
28239b2de37Sjeremylt             nnodes: total number of nodes
28339b2de37Sjeremylt             nqpts: total number of quadrature points
284d979a051Sjeremylt             *interp: Numpy array holding the row-major (nqpts * nnodes) matrix
285d979a051Sjeremylt                       expressing the values of nodal basis functions at
28639b2de37Sjeremylt                       quadrature points
287d979a051Sjeremylt             *grad: Numpy array holding the row-major (nqpts * dim * nnodes)
288d979a051Sjeremylt                     matrix expressing the derivatives of nodal basis functions
289d979a051Sjeremylt                     at quadrature points
290d979a051Sjeremylt             *qref: Numpy array of length (nqpts * dim) holding the locations of
291d979a051Sjeremylt                     quadrature points on the reference element [-1, 1]
292d979a051Sjeremylt             *qweight: Numpy array of length nnodes holding the quadrature
293d979a051Sjeremylt                        weights on the reference element
29439b2de37Sjeremylt
29539b2de37Sjeremylt           Returns:
29639b2de37Sjeremylt             basis: Ceed Basis"""
29739b2de37Sjeremylt
2987a7b0fa3SJed Brown        return BasisH1(self, topo, ncomp, nnodes, nqpts,
2997a7b0fa3SJed Brown                       interp, grad, qref, qweight)
30039b2de37Sjeremylt
30139b2de37Sjeremylt    # CeedQFunction
30239b2de37Sjeremylt    def QFunction(self, vlength, f, source):
303d979a051Sjeremylt        """Ceed QFunction: point-wise operation at quadrature points for
304d979a051Sjeremylt             evaluating volumetric terms.
30539b2de37Sjeremylt
30639b2de37Sjeremylt           Args:
30739b2de37Sjeremylt             vlength: vector length. Caller must ensure that number of quadrature
30839b2de37Sjeremylt                        points is a multiple of vlength
30939b2de37Sjeremylt             f: ctypes function pointer to evaluate action at quadrature points
31039b2de37Sjeremylt             source: absolute path to source of QFunction,
31139b2de37Sjeremylt               "\\abs_path\\file.h:function_name
31239b2de37Sjeremylt
31339b2de37Sjeremylt           Returns:
31439b2de37Sjeremylt             qfunction: Ceed QFunction"""
31539b2de37Sjeremylt
31639b2de37Sjeremylt        return QFunction(self, vlength, f, source)
31739b2de37Sjeremylt
31839b2de37Sjeremylt    def QFunctionByName(self, name):
31939b2de37Sjeremylt        """Ceed QFunction By Name: point-wise operation at quadrature points
32039b2de37Sjeremylt             from a given gallery, for evaluating volumetric terms.
32139b2de37Sjeremylt
32239b2de37Sjeremylt           Args:
32339b2de37Sjeremylt             name: name of QFunction to use from gallery
32439b2de37Sjeremylt
32539b2de37Sjeremylt           Returns:
32639b2de37Sjeremylt             qfunction: Ceed QFunction By Name"""
32739b2de37Sjeremylt
32839b2de37Sjeremylt        return QFunctionByName(self, name)
32939b2de37Sjeremylt
33039b2de37Sjeremylt    def IdentityQFunction(self, size, inmode, outmode):
33139b2de37Sjeremylt        """Ceed Idenity QFunction: identity qfunction operation.
33239b2de37Sjeremylt
33339b2de37Sjeremylt           Args:
33439b2de37Sjeremylt             size: size of the qfunction fields
33539b2de37Sjeremylt             **inmode: CeedEvalMode for input to Ceed QFunction
33639b2de37Sjeremylt             **outmode: CeedEvalMode for output to Ceed QFunction
33739b2de37Sjeremylt
33839b2de37Sjeremylt           Returns:
33939b2de37Sjeremylt             qfunction: Ceed Identity QFunction"""
34039b2de37Sjeremylt
34139b2de37Sjeremylt        return IdentityQFunction(self, size, inmode, outmode)
34239b2de37Sjeremylt
34339b2de37Sjeremylt    # CeedOperator
34439b2de37Sjeremylt    def Operator(self, qf, dqf=None, qdfT=None):
34539b2de37Sjeremylt        """Ceed Operator: composed FE-type operations on vectors.
34639b2de37Sjeremylt
34739b2de37Sjeremylt           Args:
348d979a051Sjeremylt             qf: QFunction defining the action of the operator at quadrature
349d979a051Sjeremylt                   points
35039b2de37Sjeremylt             **dqf: QFunction defining the action of the Jacobian, default None
351d979a051Sjeremylt             **dqfT: QFunction defining the action of the transpose of the
352d979a051Sjeremylt                       Jacobian, default None
35339b2de37Sjeremylt
35439b2de37Sjeremylt           Returns:
35539b2de37Sjeremylt             operator: Ceed Operator"""
35639b2de37Sjeremylt
35739b2de37Sjeremylt        return Operator(self, qf, dqf, qdfT)
35839b2de37Sjeremylt
35939b2de37Sjeremylt    def CompositeOperator(self):
36039b2de37Sjeremylt        """Composite Ceed Operator: composition of multiple CeedOperators.
36139b2de37Sjeremylt
36239b2de37Sjeremylt           Returns:
36339b2de37Sjeremylt             operator: Ceed Composite Operator"""
36439b2de37Sjeremylt
36539b2de37Sjeremylt        return CompositeOperator(self)
36639b2de37Sjeremylt
36739b2de37Sjeremylt    # Destructor
36839b2de37Sjeremylt    def __del__(self):
36939b2de37Sjeremylt        # libCEED call
37039b2de37Sjeremylt        lib.CeedDestroy(self._pointer)
37139b2de37Sjeremylt
37239b2de37Sjeremylt# ------------------------------------------------------------------------------
373