1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207
|
// This codes shows how to solve an eigenvalue problems
// in different coordinate systems. The example here is
// the Laplace equation in spherical coordinates, in the
// [theta,phi] plane. Solutions are the spherical harmonics.
// the [theta,phi] plane is a 2d domain, using SLEPc.
//
// Here, we consider the problem in spherical coordinates
// \theta \in [0,\pi] and \phi \in [0,2\pi], for a fixed
// radius: -\nabla^2 u = E u
//
// The eigenspectrum is :
// E_n = n(n+1) ; deg(n) = 2n+1
//
// Usage :
// mpirun -np 4 FreeFem++-mpi -wg laplace-2d-spherical-harmonics-SLEPc.edp \
// -split 1 -npts 400 -nev 15 -sigma 0.0
//
// Authors: Julien Garaud <julien.garaud@gmail.com>
// Pierre Jolivet <pierre.jolivet@enseeiht.fr>
/***************************************/
/* Geometry parameters */
/***************************************/
int[int] Labels = [1,2,3,4]; // labels : bottom, right, top, left sides
int[int] labPeriodic = [Labels[0],Labels[2]];
/**************************************/
/* Load PETSc & SLEPc macros */
/**************************************/
load "PETSc" // PETSc plugin
load "SLEPc" // SLEPc plugin
macro partitioner()metis// End Of Macro // metis, scotch, or parmetis
macro dimension( )2// End Of Macro // 2D or 3D
include "macro_ddm.idp" // Additional DDM functions
macro def(i)i// EOM
macro init(i)i// EOM
macro Pk() P1, periodic=[[Labels[0],x],[Labels[2],x]]//EOM
/***************************************/
/* Options for distributed solver */
/***************************************/
int s = getARGV("-split", 1) ; // Refinement factor
//
int Npts = getARGV("-npts" , 800) ; // Number of points on the perimeter
//
int nEV = getARGV("-nev" , 5) ; // Number of eigenvalues
real sigma = getARGV("-sigma", 0.0) ; // Shift
//
real radius = getARGV("-radius",1.0); // Radius of the sphere
/***************************************/
/* Verbosity and passed options */
/***************************************/
if(verbosity > 0 && mpirank == 0) {
cout << "********************************************" << endl
<< " --- " << mpirank << "/" << mpisize
<< "- laplace-2d-spherical-harmonics-SLEPc.edp " << endl
<< "********************************************" << endl
<< "- input parameters: " << endl
<< " refinement factor = " << s << endl
<< "********************************************" << endl
<< " nb of pts on perimeter = " << Npts << endl
<< "********************************************" << endl
<< " nb of eigenvalues = " << nEV << endl
<< " value of the shift = " << sigma << endl
<< "********************************************" << endl
<< " Radius of the sphere = " << radius << endl
<< "********************************************" << endl
<< endl;
}
/***************************************/
/* ############################### */
/***************************************/
meshN Th = square(1, 1); // Local mesh
int[int] arrayIntersection; // Rank of neighborings subdomains
int[int][int] restrictionIntersection(0); // Local-to-neighbors renumbering
real[int] D; // Partition of unity
/***************************************/
/* Finite Element space */
/***************************************/
// Definition of the finite element space on the domain Th
// P1 are the first order Lagrange elements
fespace Vh(Th, Pk); // local finite element space
/***************************************/
/* ############## */
/***************************************/
{ // Construction of the rectangular domain
int Thetapts = int(Npts/6.0); // pts on the x-axis sides
int Phipts = int(Npts/3.0); // pts on the y-axis sides
meshN ThBorder, ThGlobal ;
ThGlobal = square(Thetapts,Phipts,[x*pi,2.0*pi*y],label=Labels);
// .....
buildPeriodic(Th, // The local mesh
ThBorder, // The interface mesh
ThGlobal, // The global mesh
10, // Fake interface
s, // Refinement factor
1, // overlap
D, // partition of unity
arrayIntersection, // ranks of neighboring subdomains
restrictionIntersection, // local-to-neighbors renumbering
Vh, // The local Finite Element space
Pk, // FE-space
mpiCommWorld, // Communicator
false, // excluded
labPeriodic // Array of labels for periodic boundaries
);
}
/***********************************************************************/
/* Coordinate dependant differential operators */
/***********************************************************************/
/* Spherical coordinates in the (theta,phi)-plane */
/* */
/* r --> not used */
/* theta --> x in [0,pi] */
/* phi --> y in [0,2*pi] */
/* */
/* Jacobian determinant on the half-plane */
/* The det(J) = r^2*sin(theta) --> r^2*sin(x) */
/* */
macro Jac()( radius^2*sin(x) ) // End Of Macro /* */
;; /* The Jacobian */
/* */
/* The gradiant operator in spherical coordinates */
/* */
/* d/dr */
/* grad = 1/r*d/dtheta -> 1/radius*d/dx */
/* 1/(r*sin(theta))*d/dphi -> */
/* 1/(radius*sin(x)*d/dy */
/* */
macro Grad(u) [dx(u)/radius,dy(u)/(radius*sin(x))] // End Of Macro /* */
;; /* The Gradient operator */
macro Lap(u,v) ( Grad(u)'*Grad(v)) //') // End Of Macro /* */
;; /* The Laplace operator */
/* */
/* */
/***********************************************************************/
/***************************************/
/* Problem parameters */
/***************************************/
/***************************************/
/* Problem definition */
/***************************************/
varf vA(uh,vh)= intN(Th) // Definion of the problem
(Jac*(Lap(uh,vh)-sigma*uh*vh ))// Bilinear form
;
varf vB(uh,vh)= intN(Th) // Definion of the problem
( Jac*uh*vh ) // Bilinear form
;
matrix<real> A = vA(Vh,Vh);
matrix<real> B = vB(Vh,Vh);
/***************************************/
/* Build distributed matrices */
/***************************************/
dmatrix DistA(A, arrayIntersection, restrictionIntersection, D);
dmatrix DistB(DistA, B);
/***************************************/
/* Problem resolution */
/***************************************/
real[int] EigenVAL(0); // array to store eigenvalues
Vh<real>[int] def(EigenVEC)(1); // array to store eigenvectors
string ssparams = // Parameters for the distributed EigenValue solver
" -eps_nev " + nEV + // Number of eigenvalues
" -eps_type krylovschur" +
" -eps_target "+ sigma + // Shift value
" -st_type sinvert " +
" -st_pc_type lu " +
" -st_pc_factor_mat_solver_package mumps" +
" -eps_view" +
" -eps_gen_hermitian" // The problem is symmetric
;
int k = deigensolver
(DistA, // matrix OP = A − sigma*B
DistB, //
vectors = EigenVEC, // Array to store the FEM-EigenFunctions
values = EigenVAL, // Array to store the EigenValues
sparams = ssparams // Parameters for the distributed EigenValue solver
);
k=min(k,nEV); // some time the number of converged eigen value
// can be greater than nev;
/***************************************/
/* View the solution */
/***************************************/
Vh<real> Temp;
for(int i=0;i<k;i++){
if(!mpirank) cout << " Eigenvalue #"+i+" = "+EigenVAL[i]<<endl;
Temp = EigenVEC[i];
plotMPI(Th, // The local mesh
Temp[], // The local solution
"Psi("+i+") EV = "+EigenVAL[i], // Comment
Pk, // Local FE-space
def, // Macro for field definition
real, // Type
2, // 2d/3d view
1 // Wait
)
}
|