PISM, A Parallel Ice Sheet Model 2.2.1-cd005eec8 committed by Constantine Khrulev on 2025-03-07
Loading...
Searching...
No Matches
DirichletData.hh
Go to the documentation of this file.
1/* Copyright (C) 2020, 2022 PISM Authors
2 *
3 * This file is part of PISM.
4 *
5 * PISM is free software; you can redistribute it and/or modify it under the
6 * terms of the GNU General Public License as published by the Free Software
7 * Foundation; either version 3 of the License, or (at your option) any later
8 * version.
9 *
10 * PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12 * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13 * details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with PISM; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19#ifndef PISM_DIRICHLETDATA_H
20#define PISM_DIRICHLETDATA_H
21
22#include "pism/util/fem/FEM.hh"
23#include <petscmat.h>
24#include "pism/util/Vector2d.hh"
25
26namespace pism {
27
28namespace array {
29class Array;
30class Scalar;
31class Vector;
32}
33
34
35namespace fem {
36
37class Element2;
38
39//* Parts shared by scalar and 2D vector Dirichlet data classes.
41public:
42 void constrain(Element2 &element);
43 operator bool() {
44 return m_indices != NULL;
45 }
46protected:
49
50 void init(const array::Scalar *indices, const array::Array *values, double weight = 1.0);
51 void finish(const array::Array *values);
52
55 double m_weight;
56};
57
59public:
60 DirichletData_Scalar(const array::Scalar *indices, const array::Scalar *values,
61 double weight = 1.0);
63
64 void enforce(const Element2 &element, double* x_e);
65 void enforce_homogeneous(const Element2 &element, double* x_e);
66 void fix_residual(double const *const * x_global, double **r_global);
67 void fix_residual_homogeneous(double **r_global);
68 void fix_jacobian(Mat J);
69protected:
71};
72
74public:
75 DirichletData_Vector(const array::Scalar *indices, const array::Vector *values,
76 double weight);
78
79 void enforce(const Element2 &element, Vector2d* x_e);
80 void enforce_homogeneous(const Element2 &element, Vector2d* x_e);
81 void fix_residual(Vector2d const *const * x_global, Vector2d **r_global);
83 void fix_jacobian(Mat J);
84protected:
86};
87
88} // end of namespace fem
89} // end of namespace pism
90
91#endif /* PISM_DIRICHLETDATA_H */
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition Vector2d.hh:29
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition Array.hh:207
void enforce_homogeneous(const Element2 &element, double *x_e)
void fix_residual_homogeneous(double **r_global)
const array::Scalar * m_values
void fix_residual(double const *const *x_global, double **r_global)
void enforce(const Element2 &element, double *x_e)
void fix_residual(Vector2d const *const *x_global, Vector2d **r_global)
void enforce(const Element2 &element, Vector2d *x_e)
void fix_residual_homogeneous(Vector2d **r)
void enforce_homogeneous(const Element2 &element, Vector2d *x_e)
const array::Vector * m_values
double m_indices_e[q1::n_chi]
void finish(const array::Array *values)
void init(const array::Scalar *indices, const array::Array *values, double weight=1.0)
const array::Scalar * m_indices
void constrain(Element2 &element)
Constrain element, i.e. ensure that quadratures do not contribute to Dirichlet nodes by marking corre...
connected_components::details::PISMArray Array
const int n_chi
Definition FEM.hh:191