Line | Branch | Exec | Source |
---|---|---|---|
1 | ! | ||
2 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | ! SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | ! Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | ! | ||
6 | ! This file is part of SLEPc. | ||
7 | ! SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | ! | ||
10 | ! Program usage: mpiexec -n <np> ./ex20f [-n <n>] [SLEPc opts] | ||
11 | ! | ||
12 | ! Description: Simple 1-D nonlinear eigenproblem. Fortran90 equivalent of ex20.c | ||
13 | ! | ||
14 | ! The command line options are: | ||
15 | ! -n <n>, where <n> = number of grid subdivisions | ||
16 | ! | ||
17 | ! ---------------------------------------------------------------------- | ||
18 | ! Solve 1-D PDE | ||
19 | ! -u'' = lambda*u | ||
20 | ! on [0,1] subject to | ||
21 | ! u(0)=0, u'(1)=u(1)*lambda*kappa/(kappa-lambda) | ||
22 | ! ---------------------------------------------------------------------- | ||
23 | ! | ||
24 | |||
25 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
26 | ! User-defined application context | ||
27 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
28 | module UserModule | ||
29 | #include <slepc/finclude/slepcnep.h> | ||
30 | use slepcnep | ||
31 | type User | ||
32 | PetscScalar kappa | ||
33 | PetscReal h | ||
34 | end type User | ||
35 | ✗ | end module UserModule | |
36 | |||
37 | 6 | program main | |
38 | #include <slepc/finclude/slepcnep.h> | ||
39 | 6 | use UserModule | |
40 | implicit none | ||
41 | |||
42 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
43 | ! Declarations | ||
44 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
45 | ! | ||
46 | ! Variables: | ||
47 | ! nep nonlinear eigensolver context | ||
48 | ! x eigenvector | ||
49 | ! lambda eigenvalue | ||
50 | ! F,J Function and Jacobian matrices | ||
51 | ! ctx user-defined context | ||
52 | |||
53 | NEP nep | ||
54 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | Vec x, v(1) |
55 | PetscScalar lambda | ||
56 | Mat F, J | ||
57 | type(User) ctx | ||
58 | NEPType tname | ||
59 | PetscInt n, i, k, nev, its, maxit, nconv, three, one | ||
60 | PetscReal tol, norm | ||
61 | PetscScalar alpha | ||
62 | PetscMPIInt rank | ||
63 | PetscBool flg | ||
64 | PetscErrorCode ierr | ||
65 | ! Note: Any user-defined Fortran routines (such as FormJacobian) | ||
66 | ! MUST be declared as external. | ||
67 | external FormFunction, FormJacobian | ||
68 | |||
69 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
70 | ! Beginning of program | ||
71 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
72 | |||
73 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,ierr)) |
74 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) |
75 | 6 | n = 128 | |
76 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr)) |
77 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
78 | 6 | write(*,'(/A,I4)') 'Nonlinear Eigenproblem, n =',n | |
79 | endif | ||
80 | |||
81 | 6 | ctx%h = 1.0/real(n) | |
82 | 6 | ctx%kappa = 1.0 | |
83 | |||
84 | 6 | three = 3 | |
85 | 6 | one = 1 | |
86 | |||
87 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
88 | ! Create matrix data structure to hold the Function and the Jacobian | ||
89 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
90 | |||
91 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD,F,ierr)) |
92 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)) |
93 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(F,ierr)) |
94 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSeqAIJSetPreallocation(F,three,PETSC_NULL_INTEGER_ARRAY,ierr)) |
95 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatMPIAIJSetPreallocation(F,three,PETSC_NULL_INTEGER_ARRAY,one,PETSC_NULL_INTEGER_ARRAY,ierr)) |
96 | |||
97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD,J,ierr)) |
98 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)) |
99 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(J,ierr)) |
100 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSeqAIJSetPreallocation(J,three,PETSC_NULL_INTEGER_ARRAY,ierr)) |
101 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatMPIAIJSetPreallocation(J,three,PETSC_NULL_INTEGER_ARRAY,one,PETSC_NULL_INTEGER_ARRAY,ierr)) |
102 | |||
103 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
104 | ! Create the eigensolver and set various options | ||
105 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
106 | |||
107 | ! ** Create eigensolver context | ||
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPCreate(PETSC_COMM_WORLD,nep,ierr)) |
109 | |||
110 | ! ** Set routines for evaluation of Function and Jacobian | ||
111 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetFunction(nep,F,F,FormFunction,ctx,ierr)) |
112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetJacobian(nep,J,FormJacobian,ctx,ierr)) |
113 | |||
114 | ! ** Customize nonlinear solver | ||
115 | 6 | tol = 1e-9 | |
116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetTolerances(nep,tol,PETSC_CURRENT_INTEGER,ierr)) |
117 | 6 | k = 1 | |
118 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetDimensions(nep,k,PETSC_DETERMINE_INTEGER,PETSC_DETERMINE_INTEGER,ierr)) |
119 | |||
120 | ! ** Set solver parameters at runtime | ||
121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetFromOptions(nep,ierr)) |
122 | |||
123 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
124 | ! Solve the eigensystem | ||
125 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
126 | |||
127 | ! ** Evaluate initial guess | ||
128 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreateVecs(F,x,PETSC_NULL_VEC,ierr)) |
129 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecDuplicate(x,v(1),ierr)) |
130 | 6 | alpha = 1.0 | |
131 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecSet(v(1),alpha,ierr)) |
132 | 6 | k = 1 | |
133 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetInitialSpace(nep,k,v,ierr)) |
134 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecDestroy(v(1),ierr)) |
135 | |||
136 | ! ** Call the solver | ||
137 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSolve(nep,ierr)) |
138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetIterationNumber(nep,its,ierr)) |
139 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
140 | 6 | write(*,'(A,I3)') ' Number of NEP iterations =',its | |
141 | endif | ||
142 | |||
143 | ! ** Optional: Get some information from the solver and display it | ||
144 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetType(nep,tname,ierr)) |
145 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
146 | 6 | write(*,'(A,A10)') ' Solution method: ',tname | |
147 | endif | ||
148 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetDimensions(nep,nev,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr)) |
149 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
150 | 6 | write(*,'(A,I4)') ' Number of requested eigenvalues:',nev | |
151 | endif | ||
152 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetTolerances(nep,tol,maxit,ierr)) |
153 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
154 | 6 | write(*,'(A,F12.9,A,I5)') ' Stopping condition: tol=',tol,', maxit=',maxit | |
155 | endif | ||
156 | |||
157 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
158 | ! Display solution and clean up | ||
159 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
160 | |||
161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetConverged(nep,nconv,ierr)) |
162 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
163 | 6 | write(*,'(A,I2/)') ' Number of converged approximate eigenpairs:',nconv | |
164 | endif | ||
165 | |||
166 | ! ** Display eigenvalues and relative errors | ||
167 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (nconv .gt. 0) then |
168 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
169 | 6 | write(*,*) ' k ||T(k)x||' | |
170 | 6 | write(*,*) '----------------- ------------------' | |
171 | endif | ||
172 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | do i=0,nconv-1 |
173 | ! ** Get converged eigenpairs: (in this example they are always real) | ||
174 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetEigenpair(nep,i,lambda,PETSC_NULL_SCALAR,x,PETSC_NULL_VEC,ierr)) |
175 | |||
176 | ! ** Compute residual norm and error | ||
177 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPComputeError(nep,i,NEP_ERROR_RELATIVE,norm,ierr)) |
178 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
12 | if (rank .eq. 0) then |
179 | 6 | write(*,'(1P,E15.4,E18.4)') PetscRealPart(lambda), norm | |
180 | endif | ||
181 | enddo | ||
182 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank .eq. 0) then |
183 | 6 | write(*,*) | |
184 | endif | ||
185 | endif | ||
186 | |||
187 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPDestroy(nep,ierr)) |
188 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(F,ierr)) |
189 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(J,ierr)) |
190 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecDestroy(x,ierr)) |
191 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
6 | PetscCallA(SlepcFinalize(ierr)) |
192 | 1 | end | |
193 | |||
194 | ! --------------- Evaluate Function matrix T(lambda) ---------------- | ||
195 | |||
196 | 84 | subroutine FormFunction(nep,lambda,fun,B,ctx,ierr) | |
197 | use UserModule | ||
198 | implicit none | ||
199 | NEP nep | ||
200 | PetscScalar lambda, A(3), c, d | ||
201 | Mat fun,B | ||
202 | type(User) ctx | ||
203 | PetscReal h | ||
204 | PetscInt i, n, j(3), Istart, Iend, one, two, three | ||
205 | PetscErrorCode ierr | ||
206 | |||
207 | ! ** Compute Function entries and insert into matrix | ||
208 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
84 | PetscCall(MatGetSize(fun,n,PETSC_NULL_INTEGER,ierr)) |
209 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
84 | PetscCall(MatGetOwnershipRange(fun,Istart,Iend,ierr)) |
210 | 84 | h = ctx%h | |
211 | 84 | c = ctx%kappa/(lambda-ctx%kappa) | |
212 | 84 | d = n | |
213 | 84 | one = 1 | |
214 | 84 | two = 2 | |
215 | 84 | three = 3 | |
216 | |||
217 | ! ** Boundary points | ||
218 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
84 | if (Istart .eq. 0) then |
219 | 84 | i = 0 | |
220 | 84 | j(1) = 0 | |
221 | 84 | j(2) = 1 | |
222 | 84 | A(1) = 2.0*(d-lambda*h/3.0) | |
223 | 84 | A(2) = -d-lambda*h/6.0 | |
224 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
168 | PetscCall(MatSetValues(fun,one,[i],two,j,A,INSERT_VALUES,ierr)) |
225 | 84 | Istart = Istart + 1 | |
226 | endif | ||
227 | |||
228 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
84 | if (Iend .eq. n) then |
229 | 84 | i = n-1 | |
230 | 84 | j(1) = n-2 | |
231 | 84 | j(2) = n-1 | |
232 | 84 | A(1) = -d-lambda*h/6.0 | |
233 | 84 | A(2) = d-lambda*h/3.0+lambda*c | |
234 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
168 | PetscCall(MatSetValues(fun,one,[i],two,j,A,INSERT_VALUES,ierr)) |
235 | 84 | Iend = Iend - 1 | |
236 | endif | ||
237 | |||
238 | ! ** Interior grid points | ||
239 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
10668 | do i=Istart,Iend-1 |
240 | 10584 | j(1) = i-1 | |
241 | 10584 | j(2) = i | |
242 | 10584 | j(3) = i+1 | |
243 | 10584 | A(1) = -d-lambda*h/6.0 | |
244 | 10584 | A(2) = 2.0*(d-lambda*h/3.0) | |
245 | 10584 | A(3) = -d-lambda*h/6.0 | |
246 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
21252 | PetscCall(MatSetValues(fun,one,[i],three,j,A,INSERT_VALUES,ierr)) |
247 | enddo | ||
248 | |||
249 | ! ** Assemble matrix | ||
250 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
84 | PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY,ierr)) |
251 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
84 | PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY,ierr)) |
252 | |||
253 | end | ||
254 | |||
255 | ! --------------- Evaluate Jacobian matrix T'(lambda) --------------- | ||
256 | |||
257 | 66 | subroutine FormJacobian(nep,lambda,jac,ctx,ierr) | |
258 | use UserModule | ||
259 | implicit none | ||
260 | NEP nep | ||
261 | PetscScalar lambda, A(3), c | ||
262 | Mat jac | ||
263 | type(User) ctx | ||
264 | PetscReal h | ||
265 | PetscInt i, n, j(3), Istart, Iend, one, two, three | ||
266 | PetscErrorCode ierr | ||
267 | |||
268 | ! ** Compute Jacobian entries and insert into matrix | ||
269 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
66 | PetscCall(MatGetSize(jac,n,PETSC_NULL_INTEGER,ierr)) |
270 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
66 | PetscCall(MatGetOwnershipRange(jac,Istart,Iend,ierr)) |
271 | 66 | h = ctx%h | |
272 | 66 | c = ctx%kappa/(lambda-ctx%kappa) | |
273 | 66 | one = 1 | |
274 | 66 | two = 2 | |
275 | 66 | three = 3 | |
276 | |||
277 | ! ** Boundary points | ||
278 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
66 | if (Istart .eq. 0) then |
279 | 66 | i = 0 | |
280 | 66 | j(1) = 0 | |
281 | 66 | j(2) = 1 | |
282 | 66 | A(1) = -2.0*h/3.0 | |
283 | 66 | A(2) = -h/6.0 | |
284 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
132 | PetscCall(MatSetValues(jac,one,[i],two,j,A,INSERT_VALUES,ierr)) |
285 | 66 | Istart = Istart + 1 | |
286 | endif | ||
287 | |||
288 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
66 | if (Iend .eq. n) then |
289 | 66 | i = n-1 | |
290 | 66 | j(1) = n-2 | |
291 | 66 | j(2) = n-1 | |
292 | 66 | A(1) = -h/6.0 | |
293 | 66 | A(2) = -h/3.0-c*c | |
294 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
132 | PetscCall(MatSetValues(jac,one,[i],two,j,A,INSERT_VALUES,ierr)) |
295 | 66 | Iend = Iend - 1 | |
296 | endif | ||
297 | |||
298 | ! ** Interior grid points | ||
299 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
8382 | do i=Istart,Iend-1 |
300 | 8316 | j(1) = i-1 | |
301 | 8316 | j(2) = i | |
302 | 8316 | j(3) = i+1 | |
303 | 8316 | A(1) = -h/6.0 | |
304 | 8316 | A(2) = -2.0*h/3.0 | |
305 | 8316 | A(3) = -h/6.0 | |
306 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
16698 | PetscCall(MatSetValues(jac,one,[i],three,j,A,INSERT_VALUES,ierr)) |
307 | enddo | ||
308 | |||
309 | ! ** Assemble matrix | ||
310 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
66 | PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)) |
311 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
66 | PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)) |
312 | |||
313 | end | ||
314 | |||
315 | !/*TEST | ||
316 | ! | ||
317 | ! test: | ||
318 | ! suffix: 1 | ||
319 | ! args: -nep_target 4 | ||
320 | ! filter: sed -e "s/[0-9]\.[0-9]*E-[0-9]*/removed/g" -e "s/ Number of NEP iterations = [ 0-9]*/ Number of NEP iterations = /" | ||
321 | ! requires: !single | ||
322 | ! | ||
323 | !TEST*/ | ||
324 |