1 | !The Full Polymorphic Package |
---|
2 | !Copyright (C) Etienne Forest |
---|
3 | |
---|
4 | |
---|
5 | MODULE TPSA |
---|
6 | !use newda |
---|
7 | use definition |
---|
8 | use file_handler |
---|
9 | IMPLICIT NONE |
---|
10 | public |
---|
11 | integer,private::ndel ,nd2par,nd2part,nd2partt |
---|
12 | integer,private,dimension(lnv)::jfil,jfilt |
---|
13 | |
---|
14 | private equal,DAABSEQUAL,Dequaldacon ,equaldacon ,Iequaldacon !,AABSEQUAL 2002.10.17 |
---|
15 | private pow,powr,powr8,dlogt, GETORDER,CUTORDER,getchar,GETint |
---|
16 | private getdiff,getdATRA ,mul,dmulsc,dscmul |
---|
17 | private mulsc,scmul,imulsc,iscmul |
---|
18 | private div,ddivsc,dscdiv,divsc,scdiv,idivsc,iscdiv |
---|
19 | private unaryADD,add,daddsc,dscadd,addsc,scadd,iaddsc,iscadd |
---|
20 | private unarySUB,subs,dsubsc,dscsub,subsc,scsub,isubsc,iscsub |
---|
21 | private allocda,KILLda,A_OPT,K_opt |
---|
22 | private dexpt,dcost,dsint,dsqrtt,dtant,datanht,dtanht |
---|
23 | PRIVATE GETCHARnd2,GETintnd2,dputchar,dputint, filter,check_j,dsinHt,dCOSHt |
---|
24 | private GETintnd2t |
---|
25 | PRIVATE DEQUAL,REQUAL,varf,varf001 !,CHARINT |
---|
26 | ! PUBLIC VAR,ASS |
---|
27 | private pbbra,full_absT,asstaylor,getcharnd2s,GETintnd2s,GETintk |
---|
28 | private shiftda,shift000 |
---|
29 | !PRIVATE null_0,ALLOC_U,FILL_N,REFILL_N |
---|
30 | ! public, alloc_uni, null_uni, fill_uni, refill_uni |
---|
31 | |
---|
32 | private fill_uni_r ! new sagan |
---|
33 | |
---|
34 | private NO,ND,ND2,NP,NDPT,NV |
---|
35 | integer NP,NO,ND,ND2,NDPT,NV |
---|
36 | integer, TARGET :: NSPIN=0 |
---|
37 | integer, TARGET :: SPIN_pos=0 |
---|
38 | private old |
---|
39 | logical(lp) old |
---|
40 | logical(lp),target :: real_warning =.true. |
---|
41 | |
---|
42 | PRIVATE null_it,Set_Up,de_Set_Up,LINE_L,RING_L,kill_DALEVEL,dealloc_DASCRATCH,set_up_level |
---|
43 | private insert_da,append_da,GETINTegrate |
---|
44 | |
---|
45 | |
---|
46 | type(dalevel) scratchda(ndumt) !scratch levels of DA using linked list |
---|
47 | |
---|
48 | |
---|
49 | |
---|
50 | INTERFACE assignment (=) |
---|
51 | MODULE PROCEDURE EQUAL |
---|
52 | ! MODULE PROCEDURE DAABSEQUAL ! remove 2002.10.17 |
---|
53 | ! MODULE PROCEDURE AABSEQUAL ! remove 2002.10.17 |
---|
54 | MODULE PROCEDURE DEQUAL ! added 2002.10.17 ! check2002.10.17 |
---|
55 | MODULE PROCEDURE REQUAL ! added 2002.10.17 ! check2002.10.17 |
---|
56 | MODULE PROCEDURE Dequaldacon |
---|
57 | MODULE PROCEDURE equaldacon |
---|
58 | MODULE PROCEDURE Iequaldacon |
---|
59 | ! UNIVERSAL_TAYLOR |
---|
60 | |
---|
61 | MODULE PROCEDURE fill_uni_r |
---|
62 | MODULE PROCEDURE null_uni |
---|
63 | MODULE PROCEDURE fill_uni ! new sagan |
---|
64 | MODULE PROCEDURE refill_uni |
---|
65 | end INTERFACE |
---|
66 | |
---|
67 | |
---|
68 | |
---|
69 | INTERFACE print |
---|
70 | MODULE PROCEDURE printunitaylor |
---|
71 | END INTERFACE |
---|
72 | |
---|
73 | |
---|
74 | |
---|
75 | |
---|
76 | !@ <table border="4" cellspacing="1" bordercolor="#000000" id="AutoNumber2" width="400" height="135"> |
---|
77 | !@ <tr> |
---|
78 | !@ <td width="77" height="20" align="center"> |
---|
79 | !@ <span style="text-transform: uppercase"> |
---|
80 | !@ <font face="Times New Roman" size="1">+</font></span></td> |
---|
81 | !@ <td width="77" height="20" align="center"> |
---|
82 | !@ <span style="text-transform: uppercase"> |
---|
83 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
84 | !@ <td width="77" height="20" align="center"> |
---|
85 | !@ <span style="text-transform: uppercase"> |
---|
86 | !@ <font face="Times New Roman" size="1"> |
---|
87 | !@ Real(dp)</font></span></td> |
---|
88 | !@ <td width="78" height="20" align="center"> |
---|
89 | !@ <span style="text-transform: uppercase"> |
---|
90 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
91 | !@ <td width="56" height="20" align="center"> |
---|
92 | !@ <span style="text-transform: uppercase"> |
---|
93 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
94 | !@ </tr> |
---|
95 | !@ <tr> |
---|
96 | !@ <td width="77" height="20" align="center"> |
---|
97 | !@ <span style="text-transform: uppercase"> |
---|
98 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
99 | !@ <td width="77" height="20" align="center"> |
---|
100 | !@ <span style="text-transform: uppercase; font-weight:700"> |
---|
101 | !@ <font face="Times New Roman" size="1"> |
---|
102 | !@ <a href="i_tpsa.htm#ADD" style="text-decoration: none">add</a></font></span></td> |
---|
103 | !@ <td width="77" height="20" align="center"> |
---|
104 | !@ <span style="text-transform: uppercase; font-weight:700"> |
---|
105 | !@ <font face="Times New Roman" size="1"> |
---|
106 | !@ <a href="i_tpsa.htm#DADDSC" style="text-decoration: none">daddsc</a></font></span></td> |
---|
107 | !@ <td width="78" height="20" align="center"><b> |
---|
108 | !@ <font size="1" face="Times New Roman"> |
---|
109 | !@ <a href="i_tpsa.htm#ADDSC" style="text-decoration: none">ADDSC</a></font></b></td> |
---|
110 | !@ <td width="56" height="20" align="center"><b> |
---|
111 | !@ <font size="1" face="Times New Roman"> |
---|
112 | !@ <a href="i_tpsa.htm#IADDSC" style="text-decoration: none"> |
---|
113 | !@ IADDSC</a></font></b></td> |
---|
114 | !@ </tr> |
---|
115 | !@ <tr> |
---|
116 | !@ <td width="77" height="20" align="center"> |
---|
117 | !@ <span style="text-transform: uppercase"> |
---|
118 | !@ <font face="Times New Roman" size="1"> |
---|
119 | !@ Real(dp)</font></span></td> |
---|
120 | !@ <td width="77" height="20" align="center"> |
---|
121 | !@ <span style="text-transform: uppercase; font-weight:700"> |
---|
122 | !@ <font face="Times New Roman" size="1"> |
---|
123 | !@ <a href="i_tpsa.htm#DSCADD" style="text-decoration: none">dscadd</a></font></span></td> |
---|
124 | !@ <td width="77" height="20" align="center"> |
---|
125 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
126 | !@ <td width="78" height="20" align="center"> |
---|
127 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
128 | !@ <td width="56" height="20" align="center"> |
---|
129 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
130 | !@ </tr> |
---|
131 | !@ <tr> |
---|
132 | !@ <td width="77" height="20" align="center"> |
---|
133 | !@ <span style="text-transform: uppercase"> |
---|
134 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
135 | !@ <td width="77" height="20" align="center"><b> |
---|
136 | !@ <font size="1" face="Times New Roman"> |
---|
137 | !@ <a href="i_tpsa.htm#SCADD" style="text-decoration: none">SCADD</a></font></b></td> |
---|
138 | !@ <td width="77" height="20" align="center"> |
---|
139 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
140 | !@ <td width="78" height="20" align="center"> |
---|
141 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
142 | !@ <td width="56" height="20" align="center"> |
---|
143 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
144 | !@ </tr> |
---|
145 | !@ <tr> |
---|
146 | !@ <td width="77" height="20" align="center"> |
---|
147 | !@ <span style="text-transform: uppercase"> |
---|
148 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
149 | !@ <td width="77" height="20" align="center"><b> |
---|
150 | !@ <font size="1" face="Times New Roman"> |
---|
151 | !@ <a href="i_tpsa.htm#ISCADD" style="text-decoration: none"> |
---|
152 | !@ ISCADD</a></font></b></td> |
---|
153 | !@ <td width="77" height="20" align="center"> |
---|
154 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
155 | !@ <td width="78" height="20" align="center"> |
---|
156 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
157 | !@ <td width="56" height="20" align="center"> |
---|
158 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
159 | !@ </tr> |
---|
160 | !@ </table> |
---|
161 | INTERFACE OPERATOR (+) |
---|
162 | MODULE PROCEDURE unaryADD !@2 This is a unary operation |
---|
163 | MODULE PROCEDURE add |
---|
164 | MODULE PROCEDURE daddsc |
---|
165 | MODULE PROCEDURE dscadd |
---|
166 | MODULE PROCEDURE addsc |
---|
167 | MODULE PROCEDURE scadd |
---|
168 | MODULE PROCEDURE iaddsc |
---|
169 | MODULE PROCEDURE iscadd |
---|
170 | END INTERFACE |
---|
171 | |
---|
172 | |
---|
173 | |
---|
174 | |
---|
175 | !@ <table border="4" cellspacing="1" bordercolor="#000000" id="AutoNumber1" width="400" height="135"> |
---|
176 | !@ <tr> |
---|
177 | !@ <td width="77" height="20" align="center"> |
---|
178 | !@ <span style="text-transform: uppercase"> |
---|
179 | !@ <font face="Times New Roman" size="1">-</font></span></td> |
---|
180 | !@ <td width="77" height="20" align="center"> |
---|
181 | !@ <span style="text-transform: uppercase"> |
---|
182 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
183 | !@ <td width="77" height="20" align="center"> |
---|
184 | !@ <span style="text-transform: uppercase"> |
---|
185 | !@ <font face="Times New Roman" size="1"> |
---|
186 | !@ Real(dp)</font></span></td> |
---|
187 | !@ <td width="78" height="20" align="center"> |
---|
188 | !@ <span style="text-transform: uppercase"> |
---|
189 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
190 | !@ <td width="56" height="20" align="center"> |
---|
191 | !@ <span style="text-transform: uppercase"> |
---|
192 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
193 | !@ </tr> |
---|
194 | !@ <tr> |
---|
195 | !@ <td width="77" height="20" align="center"> |
---|
196 | !@ <span style="text-transform: uppercase"> |
---|
197 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
198 | !@ <td width="77" height="20" align="center"> |
---|
199 | !@ <span style="text-transform: uppercase"> |
---|
200 | !@ <font face="Times New Roman" size="1"> |
---|
201 | !@ <a href="i_tpsa.htm#SUBS" style="text-decoration: none; font-weight: 700">SUBS</a></font></span></td> |
---|
202 | !@ <td width="77" height="20" align="center"> |
---|
203 | !@ <span style="text-transform: uppercase"> |
---|
204 | !@ <font face="Times New Roman" size="1"> |
---|
205 | !@ <a href="i_tpsa.htm#DSUBSC" style="text-decoration: none; font-weight: 700">dSUBsc</a></font></span></td> |
---|
206 | !@ <td width="78" height="20" align="center"> |
---|
207 | !@ <font size="1" face="Times New Roman"> |
---|
208 | !@ <a href="i_tpsa.htm#SUBSC" style="text-decoration: none; font-weight: 700">SUBSC</a></font></td> |
---|
209 | !@ <td width="56" height="20" align="center"> |
---|
210 | !@ <font size="1" face="Times New Roman"> |
---|
211 | !@ <a href="i_tpsa.htm#ISUBSC" style="text-decoration: none; font-weight: 700">ISUBSC</a></font></td> |
---|
212 | !@ </tr> |
---|
213 | !@ <tr> |
---|
214 | !@ <td width="77" height="20" align="center"> |
---|
215 | !@ <span style="text-transform: uppercase"> |
---|
216 | !@ <font face="Times New Roman" size="1"> |
---|
217 | !@ Real(dp)</font></span></td> |
---|
218 | !@ <td width="77" height="20" align="center"> |
---|
219 | !@ <span style="text-transform: uppercase"> |
---|
220 | !@ <font face="Times New Roman" size="1"> |
---|
221 | !@ <a href="i_tpsa.htm#DSCSUB" style="text-decoration: none; font-weight: 700">dscSUB</a></font></span></td> |
---|
222 | !@ <td width="77" height="20" align="center"> |
---|
223 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
224 | !@ <td width="78" height="20" align="center"> |
---|
225 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
226 | !@ <td width="56" height="20" align="center"> |
---|
227 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
228 | !@ </tr> |
---|
229 | !@ <tr> |
---|
230 | !@ <td width="77" height="20" align="center"> |
---|
231 | !@ <span style="text-transform: uppercase"> |
---|
232 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
233 | !@ <td width="77" height="20" align="center"> |
---|
234 | !@ <font size="1" face="Times New Roman"> |
---|
235 | !@ <a href="i_tpsa.htm#SCSUB" style="text-decoration: none; font-weight: 700"> |
---|
236 | !@ SCSUB</a></font></td> |
---|
237 | !@ <td width="77" height="20" align="center"> |
---|
238 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
239 | !@ <td width="78" height="20" align="center"> |
---|
240 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
241 | !@ <td width="56" height="20" align="center"> |
---|
242 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
243 | !@ </tr> |
---|
244 | !@ <tr> |
---|
245 | !@ <td width="77" height="20" align="center"> |
---|
246 | !@ <span style="text-transform: uppercase"> |
---|
247 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
248 | !@ <td width="77" height="20" align="center"> |
---|
249 | !@ <font size="1" face="Times New Roman"> |
---|
250 | !@ <a href="i_tpsa.htm#ISCSUB" style="text-decoration: none; font-weight: 700"> |
---|
251 | !@ ISCSUB</a></font></td> |
---|
252 | !@ <td width="77" height="20" align="center"> |
---|
253 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
254 | !@ <td width="78" height="20" align="center"> |
---|
255 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
256 | !@ <td width="56" height="20" align="center"> |
---|
257 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
258 | !@ </tr> |
---|
259 | !@ </table> |
---|
260 | INTERFACE OPERATOR (-) |
---|
261 | MODULE PROCEDURE unarySUB |
---|
262 | MODULE PROCEDURE subs |
---|
263 | MODULE PROCEDURE dsubsc |
---|
264 | MODULE PROCEDURE dscsub |
---|
265 | MODULE PROCEDURE subsc |
---|
266 | MODULE PROCEDURE scsub |
---|
267 | MODULE PROCEDURE isubsc |
---|
268 | MODULE PROCEDURE iscsub |
---|
269 | END INTERFACE |
---|
270 | |
---|
271 | |
---|
272 | |
---|
273 | !@ <table border="4" cellspacing="1" bordercolor="#000000" id="AutoNumber1" width="400" height="134"> |
---|
274 | !@ <tr> |
---|
275 | !@ <td width="77" height="20" align="center"> |
---|
276 | !@ <span style="text-transform: uppercase"> |
---|
277 | !@ <font face="Times New Roman" size="1">*</font></span></td> |
---|
278 | !@ <td width="77" height="20" align="center"> |
---|
279 | !@ <span style="text-transform: uppercase"> |
---|
280 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
281 | !@ <td width="77" height="20" align="center"> |
---|
282 | !@ <span style="text-transform: uppercase"> |
---|
283 | !@ <font face="Times New Roman" size="1"> |
---|
284 | !@ Real(dp)</font></span></td> |
---|
285 | !@ <td width="78" height="20" align="center"> |
---|
286 | !@ <span style="text-transform: uppercase"> |
---|
287 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
288 | !@ <td width="56" height="20" align="center"> |
---|
289 | !@ <span style="text-transform: uppercase"> |
---|
290 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
291 | !@ </tr> |
---|
292 | !@ <tr> |
---|
293 | !@ <td width="77" height="20" align="center"> |
---|
294 | !@ <span style="text-transform: uppercase"> |
---|
295 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
296 | !@ <td width="77" height="20" align="center"> |
---|
297 | !@ <span style="text-transform: uppercase"> |
---|
298 | !@ <font face="Times New Roman" size="1"> |
---|
299 | !@ <a href="i_tpsa.htm#MUL" style="text-decoration: none; font-weight:700">MUL</a></font></span></td> |
---|
300 | !@ <td width="77" height="20" align="center"> |
---|
301 | !@ <span style="text-transform: uppercase"> |
---|
302 | !@ <font face="Times New Roman" size="1"> |
---|
303 | !@ <a href="i_tpsa.htm#DMULSC" style="text-decoration: none; font-weight:700">dMULsc</a></font></span></td> |
---|
304 | !@ <td width="78" height="20" align="center"> |
---|
305 | !@ <font size="1" face="Times New Roman"> |
---|
306 | !@ <a href="i_tpsa.htm#MULSC" style="text-decoration: none; font-weight:700">MULSC</a></font></td> |
---|
307 | !@ <td width="56" height="20" align="center"> |
---|
308 | !@ <font size="1" face="Times New Roman"> |
---|
309 | !@ <a href="i_tpsa.htm#IMULSC" style="text-decoration: none; font-weight:700">IMULSC</a></font></td> |
---|
310 | !@ </tr> |
---|
311 | !@ <tr> |
---|
312 | !@ <td width="77" height="19" align="center"> |
---|
313 | !@ <span style="text-transform: uppercase"> |
---|
314 | !@ <font face="Times New Roman" size="1"> |
---|
315 | !@ Real(dp)</font></span></td> |
---|
316 | !@ <td width="77" height="19" align="center"> |
---|
317 | !@ <span style="text-transform: uppercase"> |
---|
318 | !@ <font face="Times New Roman" size="1"> |
---|
319 | !@ <a href="i_tpsa.htm#DSCMUL" style="text-decoration: none; font-weight:700">dscMUL</a></font></span></td> |
---|
320 | !@ <td width="77" height="19" align="center"> |
---|
321 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
322 | !@ <td width="78" height="19" align="center"> |
---|
323 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
324 | !@ <td width="56" height="19" align="center"> |
---|
325 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
326 | !@ </tr> |
---|
327 | !@ <tr> |
---|
328 | !@ <td width="77" height="20" align="center"> |
---|
329 | !@ <span style="text-transform: uppercase"> |
---|
330 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
331 | !@ <td width="77" height="20" align="center"> |
---|
332 | !@ <font size="1" face="Times New Roman"> |
---|
333 | !@ <a href="i_tpsa.htm#SCMUL" style="text-decoration: none; font-weight:700"> |
---|
334 | !@ SCMUL</a></font></td> |
---|
335 | !@ <td width="77" height="20" align="center"> |
---|
336 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
337 | !@ <td width="78" height="20" align="center"> |
---|
338 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
339 | !@ <td width="56" height="20" align="center"> |
---|
340 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
341 | !@ </tr> |
---|
342 | !@ <tr> |
---|
343 | !@ <td width="77" height="20" align="center"> |
---|
344 | !@ <span style="text-transform: uppercase"> |
---|
345 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
346 | !@ <td width="77" height="20" align="center"> |
---|
347 | !@ <font size="1" face="Times New Roman"> |
---|
348 | !@ <a href="i_tpsa.htm#ISCMUL" style="text-decoration: none; font-weight:700"> |
---|
349 | !@ ISCMUL</a></font></td> |
---|
350 | !@ <td width="77" height="20" align="center"> |
---|
351 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
352 | !@ <td width="78" height="20" align="center"> |
---|
353 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
354 | !@ <td width="56" height="20" align="center"> |
---|
355 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
356 | !@ </tr> |
---|
357 | !@ </table> |
---|
358 | |
---|
359 | INTERFACE OPERATOR (*) |
---|
360 | MODULE PROCEDURE mul |
---|
361 | MODULE PROCEDURE dmulsc |
---|
362 | MODULE PROCEDURE dscmul |
---|
363 | MODULE PROCEDURE mulsc |
---|
364 | MODULE PROCEDURE scmul |
---|
365 | MODULE PROCEDURE imulsc |
---|
366 | MODULE PROCEDURE iscmul |
---|
367 | END INTERFACE |
---|
368 | |
---|
369 | !@ <table border="4" cellspacing="1" bordercolor="#000000" id="AutoNumber1" width="400" height="135"> |
---|
370 | !@ <tr> |
---|
371 | !@ <td width="0" height="20" align="center"> |
---|
372 | !@ <span style="text-transform: uppercase"> |
---|
373 | !@ <font face="Times New Roman" size="1">/</font></span></td> |
---|
374 | !@ <td width="0" height="20" align="center"> |
---|
375 | !@ <span style="text-transform: uppercase"> |
---|
376 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
377 | !@ <td width="0" height="20" align="center"> |
---|
378 | !@ <span style="text-transform: uppercase"> |
---|
379 | !@ <font face="Times New Roman" size="1"> |
---|
380 | !@ Real(dp)</font></span></td> |
---|
381 | !@ <td width="0" height="20" align="center"> |
---|
382 | !@ <span style="text-transform: uppercase"> |
---|
383 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
384 | !@ <td width="0" height="20" align="center"> |
---|
385 | !@ <span style="text-transform: uppercase"> |
---|
386 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
387 | !@ </tr> |
---|
388 | !@ <tr> |
---|
389 | !@ <td width="0" height="20" align="center"> |
---|
390 | !@ <span style="text-transform: uppercase"> |
---|
391 | !@ <font face="Times New Roman" size="1">Taylor</font></span></td> |
---|
392 | !@ <td width="0" height="20" align="center"> |
---|
393 | !@ <span style="text-transform: uppercase"> |
---|
394 | !@ <font face="Times New Roman" size="1"> |
---|
395 | !@ <a href="i_tpsa.htm#DIV" style="text-decoration: none; font-weight: 700">div</a></font></span></td> |
---|
396 | !@ <td width="0" height="20" align="center"> |
---|
397 | !@ <span style="text-transform: uppercase"> |
---|
398 | !@ <font face="Times New Roman" size="1"> |
---|
399 | !@ <a href="i_tpsa.htm#DDIVSC" style="text-decoration: none; font-weight: 700">dDIVsc</a></font></span></td> |
---|
400 | !@ <td width="0" height="20" align="center"><font size="1"> |
---|
401 | !@ <a href="i_tpsa.htm#DIVSC" style="text-decoration: none; font-weight: 700">DIVSC</a></font></td> |
---|
402 | !@ <td width="0" height="20" align="center"><font size="1"> |
---|
403 | !@ <a href="i_tpsa.htm#IDIVSC" style="text-decoration: none; font-weight: 700"> |
---|
404 | !@ IDIVSC</a></font></td> |
---|
405 | !@ </tr> |
---|
406 | !@ <tr> |
---|
407 | !@ <td width="0" height="20" align="center"> |
---|
408 | !@ <span style="text-transform: uppercase"> |
---|
409 | !@ <font face="Times New Roman" size="1"> |
---|
410 | !@ Real(dp)</font></span></td> |
---|
411 | !@ <td width="0" height="20" align="center"> |
---|
412 | !@ <span style="text-transform: uppercase"> |
---|
413 | !@ <font face="Times New Roman" size="1"> |
---|
414 | !@ <a href="i_tpsa.htm#DSCDIV" style="text-decoration: none; font-weight: 700">dscDIV</a></font></span></td> |
---|
415 | !@ <td width="0" height="20" align="center"> |
---|
416 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
417 | !@ <td width="0" height="20" align="center"> |
---|
418 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
419 | !@ <td width="0" height="20" align="center"> |
---|
420 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
421 | !@ </tr> |
---|
422 | !@ <tr> |
---|
423 | !@ <td width="0" height="20" align="center"> |
---|
424 | !@ <span style="text-transform: uppercase"> |
---|
425 | !@ <font face="Times New Roman" size="1">Real(sp)</font></span></td> |
---|
426 | !@ <td width="0" height="20" align="center"><font size="1"> |
---|
427 | !@ <a href="i_tpsa.htm#SCDIV" style="text-decoration: none; font-weight: 700"> |
---|
428 | !@ SCDIV</a></font></td> |
---|
429 | !@ <td width="0" height="20" align="center"> |
---|
430 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
431 | !@ <td width="0" height="20" align="center"> |
---|
432 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
433 | !@ <td width="0" height="20" align="center"> |
---|
434 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
435 | !@ </tr> |
---|
436 | !@ <tr> |
---|
437 | !@ <td width="0" height="20" align="center"> |
---|
438 | !@ <span style="text-transform: uppercase"> |
---|
439 | !@ <font face="Times New Roman" size="1">Integer</font></span></td> |
---|
440 | !@ <td width="0" height="20" align="center"><font size="1"> |
---|
441 | !@ <a href="i_tpsa.htm#ISCDIV" style="text-decoration: none; font-weight: 700"> |
---|
442 | !@ ISCDIV</a></font></td> |
---|
443 | !@ <td width="0" height="20" align="center"> |
---|
444 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
445 | !@ <td width="0" height="20" align="center"> |
---|
446 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
447 | !@ <td width="0" height="20" align="center"> |
---|
448 | !@ <font size="2" face="Times New Roman"><b>F90</b></font></td> |
---|
449 | !@ </tr> |
---|
450 | !@ </table> |
---|
451 | |
---|
452 | INTERFACE OPERATOR (/) |
---|
453 | MODULE PROCEDURE div |
---|
454 | MODULE PROCEDURE ddivsc |
---|
455 | MODULE PROCEDURE dscdiv |
---|
456 | MODULE PROCEDURE divsc |
---|
457 | MODULE PROCEDURE scdiv |
---|
458 | MODULE PROCEDURE idivsc |
---|
459 | MODULE PROCEDURE iscdiv |
---|
460 | END INTERFACE |
---|
461 | |
---|
462 | |
---|
463 | INTERFACE OPERATOR (**) |
---|
464 | MODULE PROCEDURE POW |
---|
465 | MODULE PROCEDURE POWR |
---|
466 | MODULE PROCEDURE POWR8 |
---|
467 | END INTERFACE |
---|
468 | |
---|
469 | |
---|
470 | ! New Operators |
---|
471 | |
---|
472 | INTERFACE OPERATOR (.mono.) |
---|
473 | MODULE PROCEDURE dputint0 !@1 single integer |
---|
474 | MODULE PROCEDURE dputint !@1 Accepts J(nv) |
---|
475 | MODULE PROCEDURE dputchar !@1 Accepts String such as '12' |
---|
476 | END INTERFACE |
---|
477 | |
---|
478 | INTERFACE OPERATOR (.var.) |
---|
479 | MODULE PROCEDURE varf !@1 replaces var (overloads DAVAR) |
---|
480 | MODULE PROCEDURE varf001 !@1 replaces var001 |
---|
481 | END INTERFACE |
---|
482 | |
---|
483 | INTERFACE OPERATOR (.d.) |
---|
484 | MODULE PROCEDURE getdiff !@1 takes derivatives |
---|
485 | END INTERFACE |
---|
486 | |
---|
487 | INTERFACE OPERATOR (.i.) |
---|
488 | MODULE PROCEDURE GETINTegrate !@1 takes derivatives |
---|
489 | END INTERFACE |
---|
490 | |
---|
491 | |
---|
492 | INTERFACE OPERATOR (.SUB.) |
---|
493 | MODULE PROCEDURE GETORDER |
---|
494 | MODULE PROCEDURE getchar |
---|
495 | MODULE PROCEDURE GETint |
---|
496 | END INTERFACE |
---|
497 | |
---|
498 | INTERFACE OPERATOR (.PAR.) |
---|
499 | MODULE PROCEDURE getcharnd2 |
---|
500 | MODULE PROCEDURE GETintnd2 |
---|
501 | END INTERFACE |
---|
502 | |
---|
503 | INTERFACE OPERATOR (.part.) |
---|
504 | MODULE PROCEDURE GETintnd2t |
---|
505 | END INTERFACE |
---|
506 | |
---|
507 | |
---|
508 | INTERFACE OPERATOR (<=) |
---|
509 | MODULE PROCEDURE getcharnd2s |
---|
510 | MODULE PROCEDURE GETintnd2s |
---|
511 | MODULE PROCEDURE GETintk |
---|
512 | END INTERFACE |
---|
513 | |
---|
514 | INTERFACE OPERATOR (.CUT.) |
---|
515 | MODULE PROCEDURE CUTORDER |
---|
516 | END INTERFACE |
---|
517 | |
---|
518 | INTERFACE OPERATOR (.K.) |
---|
519 | MODULE PROCEDURE getdATRA ! Used internally primarily |
---|
520 | END INTERFACE |
---|
521 | |
---|
522 | INTERFACE OPERATOR (.pb.) |
---|
523 | MODULE PROCEDURE pbbra |
---|
524 | END INTERFACE |
---|
525 | |
---|
526 | ! intrisic functions overloaded |
---|
527 | |
---|
528 | INTERFACE abs |
---|
529 | MODULE PROCEDURE DAABSEQUAL ! remove 2002.10.17 |
---|
530 | END INTERFACE |
---|
531 | INTERFACE dabs |
---|
532 | MODULE PROCEDURE DAABSEQUAL ! remove 2002.10.17 |
---|
533 | END INTERFACE |
---|
534 | |
---|
535 | INTERFACE exp |
---|
536 | MODULE PROCEDURE dexpt |
---|
537 | END INTERFACE |
---|
538 | INTERFACE dexp |
---|
539 | MODULE PROCEDURE dexpt |
---|
540 | END INTERFACE |
---|
541 | INTERFACE cexp |
---|
542 | MODULE PROCEDURE dexpt |
---|
543 | END INTERFACE |
---|
544 | INTERFACE cdexp |
---|
545 | MODULE PROCEDURE dexpt |
---|
546 | END INTERFACE |
---|
547 | |
---|
548 | INTERFACE cos |
---|
549 | MODULE PROCEDURE dcost |
---|
550 | END INTERFACE |
---|
551 | INTERFACE cdcos |
---|
552 | MODULE PROCEDURE dcost |
---|
553 | END INTERFACE |
---|
554 | INTERFACE dcos |
---|
555 | MODULE PROCEDURE dcost |
---|
556 | END INTERFACE |
---|
557 | INTERFACE ccos |
---|
558 | MODULE PROCEDURE dcost |
---|
559 | END INTERFACE |
---|
560 | |
---|
561 | INTERFACE cosH |
---|
562 | MODULE PROCEDURE dcosHt |
---|
563 | END INTERFACE |
---|
564 | INTERFACE dcosH |
---|
565 | MODULE PROCEDURE dcosHt |
---|
566 | END INTERFACE |
---|
567 | |
---|
568 | INTERFACE sin |
---|
569 | MODULE PROCEDURE dsint |
---|
570 | END INTERFACE |
---|
571 | INTERFACE cdsin |
---|
572 | MODULE PROCEDURE dsint |
---|
573 | END INTERFACE |
---|
574 | INTERFACE ccsin |
---|
575 | MODULE PROCEDURE dsint |
---|
576 | END INTERFACE |
---|
577 | INTERFACE dsin |
---|
578 | MODULE PROCEDURE dsint |
---|
579 | END INTERFACE |
---|
580 | |
---|
581 | INTERFACE sinH |
---|
582 | MODULE PROCEDURE dsinHt |
---|
583 | END INTERFACE |
---|
584 | INTERFACE dsinH |
---|
585 | MODULE PROCEDURE dsinHt |
---|
586 | END INTERFACE |
---|
587 | |
---|
588 | INTERFACE log |
---|
589 | MODULE PROCEDURE dlogt |
---|
590 | END INTERFACE |
---|
591 | INTERFACE dlog |
---|
592 | MODULE PROCEDURE dlogt |
---|
593 | END INTERFACE |
---|
594 | INTERFACE cdlog |
---|
595 | MODULE PROCEDURE dlogt |
---|
596 | END INTERFACE |
---|
597 | INTERFACE clog |
---|
598 | MODULE PROCEDURE dlogt |
---|
599 | END INTERFACE |
---|
600 | |
---|
601 | INTERFACE sqrt |
---|
602 | MODULE PROCEDURE dsqrtt |
---|
603 | END INTERFACE |
---|
604 | INTERFACE dsqrt |
---|
605 | MODULE PROCEDURE dsqrtt |
---|
606 | END INTERFACE |
---|
607 | |
---|
608 | INTERFACE atanh |
---|
609 | MODULE PROCEDURE datanht |
---|
610 | END INTERFACE |
---|
611 | |
---|
612 | |
---|
613 | |
---|
614 | INTERFACE tan |
---|
615 | MODULE PROCEDURE dtant |
---|
616 | END INTERFACE |
---|
617 | |
---|
618 | INTERFACE tanh |
---|
619 | MODULE PROCEDURE dtanht |
---|
620 | END INTERFACE |
---|
621 | |
---|
622 | |
---|
623 | INTERFACE dtan |
---|
624 | MODULE PROCEDURE dtant |
---|
625 | END INTERFACE |
---|
626 | |
---|
627 | ! Non-intrisic Functions |
---|
628 | |
---|
629 | INTERFACE pek |
---|
630 | MODULE PROCEDURE pek000 ! not private |
---|
631 | END INTERFACE |
---|
632 | |
---|
633 | INTERFACE pok |
---|
634 | MODULE PROCEDURE pok000 ! not private |
---|
635 | END INTERFACE |
---|
636 | |
---|
637 | INTERFACE shiftda |
---|
638 | MODULE PROCEDURE shift000 ! not private |
---|
639 | END INTERFACE |
---|
640 | |
---|
641 | ! INTERFACE var |
---|
642 | ! MODULE PROCEDURE var000 ! not private |
---|
643 | ! MODULE PROCEDURE var001 ! not private |
---|
644 | ! END INTERFACE |
---|
645 | |
---|
646 | INTERFACE cfu |
---|
647 | MODULE PROCEDURE cfu000 ! not private |
---|
648 | END INTERFACE |
---|
649 | |
---|
650 | INTERFACE full_abs |
---|
651 | MODULE PROCEDURE full_absT |
---|
652 | END INTERFACE |
---|
653 | |
---|
654 | ! INTERFACE daread |
---|
655 | ! MODULE PROCEDURE rea |
---|
656 | ! END INTERFACE |
---|
657 | |
---|
658 | ! INTERFACE read |
---|
659 | ! MODULE PROCEDURE rea |
---|
660 | ! END INTERFACE |
---|
661 | |
---|
662 | ! INTERFACE daprint |
---|
663 | ! MODULE PROCEDURE pri |
---|
664 | ! END INTERFACE |
---|
665 | |
---|
666 | |
---|
667 | |
---|
668 | ! INTERFACE print |
---|
669 | ! MODULE PROCEDURE pri |
---|
670 | ! END INTERFACE |
---|
671 | |
---|
672 | |
---|
673 | ! Constructors and Destructors |
---|
674 | |
---|
675 | INTERFACE alloc |
---|
676 | MODULE PROCEDURE allocda |
---|
677 | MODULE PROCEDURE A_OPT |
---|
678 | MODULE PROCEDURE allocdas |
---|
679 | MODULE PROCEDURE alloc_u |
---|
680 | END INTERFACE |
---|
681 | |
---|
682 | INTERFACE KILL |
---|
683 | MODULE PROCEDURE KILLda |
---|
684 | MODULE PROCEDURE KILLdas |
---|
685 | MODULE PROCEDURE K_opt |
---|
686 | MODULE PROCEDURE kill_uni |
---|
687 | END INTERFACE |
---|
688 | |
---|
689 | INTERFACE alloctpsa |
---|
690 | MODULE PROCEDURE allocda |
---|
691 | END INTERFACE |
---|
692 | |
---|
693 | INTERFACE KILLtpsa |
---|
694 | MODULE PROCEDURE KILLda |
---|
695 | END INTERFACE |
---|
696 | |
---|
697 | |
---|
698 | ! management routines |
---|
699 | |
---|
700 | INTERFACE ass |
---|
701 | MODULE PROCEDURE asstaylor !2000.12.25 |
---|
702 | END INTERFACE |
---|
703 | |
---|
704 | |
---|
705 | |
---|
706 | CONTAINS |
---|
707 | |
---|
708 | subroutine fliptaylor(xy,xyf,i) |
---|
709 | implicit none |
---|
710 | type(taylor), intent(inout) :: xy,xyf |
---|
711 | integer i |
---|
712 | call flip_i(xy%i,xyf%i,i) |
---|
713 | end subroutine fliptaylor |
---|
714 | |
---|
715 | |
---|
716 | SUBROUTINE change_default_tpsa(i) |
---|
717 | implicit none |
---|
718 | INTEGER, intent(in) :: I |
---|
719 | if(last_tpsa==0) then |
---|
720 | if(i==1) then |
---|
721 | default_tpsa=.true. |
---|
722 | if(i==1.and.lingyun_yang )write(6,*) " Default TPSA is CPP package of Yang" |
---|
723 | call change_package(i) |
---|
724 | else |
---|
725 | default_tpsa=.false. |
---|
726 | call change_package(i) |
---|
727 | if(i==2.and.(.not.lingyun_yang) )write(6,*) " Default TPSA is FORTRAN package of Berz (LBNL)" |
---|
728 | endif |
---|
729 | else |
---|
730 | write(6,*) " You could not change default TPSA here " |
---|
731 | write(6,*) " Only prior to any call to TPSA or PTC or after a PTC_END " |
---|
732 | stop 666 |
---|
733 | endif |
---|
734 | end SUBROUTINE change_default_tpsa |
---|
735 | |
---|
736 | |
---|
737 | subroutine set_in_tpsa(NO1,ND1,ND21,NP1,NDPT1,NV1,log) |
---|
738 | implicit none |
---|
739 | integer NO1,ND1,ND21,NP1,NDPT1,NV1 |
---|
740 | logical(lp) log |
---|
741 | old=log |
---|
742 | NO=NO1 |
---|
743 | ND=ND1 |
---|
744 | ND2=ND21 |
---|
745 | NP=NP1 |
---|
746 | NDPT=NDPT1 |
---|
747 | NV=NV1 |
---|
748 | end subroutine set_in_tpsa |
---|
749 | |
---|
750 | subroutine count_taylor(n,ns,ne) |
---|
751 | implicit none |
---|
752 | integer n,ns,ne,i |
---|
753 | call count_da(n) |
---|
754 | ns=0 |
---|
755 | do i=1,ndumt |
---|
756 | ns=scratchda(i)%n+ns |
---|
757 | enddo |
---|
758 | ne=n-ns |
---|
759 | end subroutine count_taylor |
---|
760 | |
---|
761 | |
---|
762 | FUNCTION unaryADD( S1 ) |
---|
763 | implicit none |
---|
764 | TYPE (TAYLOR) unaryADD |
---|
765 | TYPE (TAYLOR), INTENT (IN) :: S1 |
---|
766 | integer localmaster |
---|
767 | |
---|
768 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
769 | |
---|
770 | localmaster=master |
---|
771 | |
---|
772 | ! call check(s1) |
---|
773 | call ass(unaryADD) |
---|
774 | |
---|
775 | unaryADD=s1 |
---|
776 | |
---|
777 | master=localmaster |
---|
778 | |
---|
779 | END FUNCTION unaryADD |
---|
780 | |
---|
781 | FUNCTION unarySUB( S1 ) |
---|
782 | implicit none |
---|
783 | TYPE (TAYLOR) unarySUB |
---|
784 | TYPE (TAYLOR), INTENT (IN) :: S1 |
---|
785 | integer localmaster |
---|
786 | |
---|
787 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
788 | localmaster=master |
---|
789 | |
---|
790 | ! call check(s1) |
---|
791 | call ass(unarySUB) |
---|
792 | |
---|
793 | ! unarySUB=(-one)*s1 |
---|
794 | ! if(old) then |
---|
795 | call dacmu(s1%i,-1.0_dp,temp) |
---|
796 | call dacop(temp,unarySUB%i) |
---|
797 | ! else |
---|
798 | ! call newdacmu(s1%j,-one,unarySUB%j) |
---|
799 | ! ! call newdacmu(s1%j,-one,templ) |
---|
800 | ! ! call newdacop(templ,unarySUB%j) |
---|
801 | ! endif |
---|
802 | master=localmaster |
---|
803 | |
---|
804 | END FUNCTION unarySUB |
---|
805 | |
---|
806 | SUBROUTINE maketree(S1,s2) |
---|
807 | implicit none |
---|
808 | type (TAYLOR),INTENT(IN)::S1 |
---|
809 | type (TAYLOR),INTENT(inOUT):: s2 |
---|
810 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
811 | |
---|
812 | ! if(old) then |
---|
813 | call mtree((/s1%i/),1,(/s2%i/),1) |
---|
814 | ! else |
---|
815 | ! call newdacop(s1%j,s2%j) |
---|
816 | ! endif |
---|
817 | END SUBROUTINE maketree |
---|
818 | |
---|
819 | SUBROUTINE allocda(S1) |
---|
820 | implicit none |
---|
821 | type (TAYLOR),INTENT(INOUT)::S1 |
---|
822 | |
---|
823 | ! IF(first_time) THEN |
---|
824 | IF(last_tpsa==0) THEN |
---|
825 | w_p=0 |
---|
826 | w_p%nc=1 |
---|
827 | w_p=(/" No TPSA package ever initialized "/) |
---|
828 | w_p%fc='(1((1X,A72),/))' |
---|
829 | ! call !write_e(111) |
---|
830 | ENDIF |
---|
831 | ! if(old) then |
---|
832 | s1%i=0 |
---|
833 | call etall1(s1%i) |
---|
834 | ! else |
---|
835 | ! call nullnewda(s1%j) |
---|
836 | ! call allocnewda(s1%j) |
---|
837 | ! endif |
---|
838 | END SUBROUTINE allocda |
---|
839 | |
---|
840 | SUBROUTINE A_OPT(S1,S2,s3,s4,s5,s6,s7,s8,s9,s10) |
---|
841 | implicit none |
---|
842 | type (taylor),INTENT(INout)::S1,S2 |
---|
843 | type (taylor),optional, INTENT(INout):: s3,s4,s5,s6,s7,s8,s9,s10 |
---|
844 | call allocda(s1) |
---|
845 | call allocda(s2) |
---|
846 | if(present(s3)) call allocda(s3) |
---|
847 | if(present(s4)) call allocda(s4) |
---|
848 | if(present(s5)) call allocda(s5) |
---|
849 | if(present(s6)) call allocda(s6) |
---|
850 | if(present(s7)) call allocda(s7) |
---|
851 | if(present(s8)) call allocda(s8) |
---|
852 | if(present(s9)) call allocda(s9) |
---|
853 | if(present(s10))call allocda(s10) |
---|
854 | END SUBROUTINE A_opt |
---|
855 | |
---|
856 | SUBROUTINE K_OPT(S1,S2,s3,s4,s5,s6,s7,s8,s9,s10) |
---|
857 | implicit none |
---|
858 | type (taylor),INTENT(INout)::S1,S2 |
---|
859 | type (taylor),optional, INTENT(INout):: s3,s4,s5,s6,s7,s8,s9,s10 |
---|
860 | call KILLDA(s1) |
---|
861 | call KILLDA(s2) |
---|
862 | if(present(s3)) call KILLDA(s3) |
---|
863 | if(present(s4)) call KILLDA(s4) |
---|
864 | if(present(s5)) call KILLDA(s5) |
---|
865 | if(present(s6)) call KILLDA(s6) |
---|
866 | if(present(s7)) call KILLDA(s7) |
---|
867 | if(present(s8)) call KILLDA(s8) |
---|
868 | if(present(s9)) call KILLDA(s9) |
---|
869 | if(present(s10))call KILLDA(s10) |
---|
870 | END SUBROUTINE K_opt |
---|
871 | |
---|
872 | |
---|
873 | SUBROUTINE ALLOCDAS(S1,k) |
---|
874 | implicit none |
---|
875 | type (TAYLOR),INTENT(INOUT),dimension(:)::S1 |
---|
876 | INTEGER,optional,INTENT(IN)::k |
---|
877 | INTEGER J,i,N |
---|
878 | |
---|
879 | if(present(k)) then |
---|
880 | I=LBOUND(S1,DIM=1) |
---|
881 | N=LBOUND(S1,DIM=1)+K-1 |
---|
882 | else |
---|
883 | I=LBOUND(S1,DIM=1) |
---|
884 | N=UBOUND(S1,DIM=1) |
---|
885 | endif |
---|
886 | |
---|
887 | DO J=I,N |
---|
888 | CALL allocDA(S1(j)) |
---|
889 | ENDDO |
---|
890 | |
---|
891 | END SUBROUTINE ALLOCDAS |
---|
892 | |
---|
893 | SUBROUTINE KILLda(S1) |
---|
894 | implicit none |
---|
895 | type (TAYLOR),INTENT(INOUT)::S1 |
---|
896 | ! if(old) then |
---|
897 | call DADAL1(s1%i) |
---|
898 | ! else |
---|
899 | ! call KILLNEWDAs(s1%j) |
---|
900 | ! endif |
---|
901 | |
---|
902 | END SUBROUTINE KILLda |
---|
903 | |
---|
904 | SUBROUTINE KILLDAS(S1,k) |
---|
905 | implicit none |
---|
906 | type (TAYLOR),INTENT(INOUT),dimension(:)::S1 |
---|
907 | INTEGER,optional,INTENT(IN)::k |
---|
908 | INTEGER J,i,N |
---|
909 | |
---|
910 | if(present(k)) then |
---|
911 | I=LBOUND(S1,DIM=1) |
---|
912 | N=LBOUND(S1,DIM=1)+K-1 |
---|
913 | else |
---|
914 | I=LBOUND(S1,DIM=1) |
---|
915 | N=UBOUND(S1,DIM=1) |
---|
916 | endif |
---|
917 | |
---|
918 | DO J=I,N |
---|
919 | CALL KILLDA(S1(j)) |
---|
920 | ENDDO |
---|
921 | |
---|
922 | END SUBROUTINE KILLDAS |
---|
923 | |
---|
924 | |
---|
925 | SUBROUTINE EQUAL(S2,S1) |
---|
926 | implicit none |
---|
927 | type (TAYLOR),INTENT(inOUT)::S2 |
---|
928 | type (TAYLOR),INTENT(IN)::S1 |
---|
929 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
930 | |
---|
931 | call check_snake |
---|
932 | ! if(old) then |
---|
933 | if(s2%i==0) then |
---|
934 | call crap1("EQUAL 1 in tpsa") !call allocw(s2) |
---|
935 | endif |
---|
936 | if(s1%i==0) call crap1("EQUAL 2") ! call allocw(s1) |
---|
937 | CALL DACOP(S1%I,S2%I) |
---|
938 | ! else |
---|
939 | ! IF (.NOT. ASSOCIATED(s2%j%r)) call crap1("EQUAL 3") !call allocw(s2) |
---|
940 | ! IF (.NOT. ASSOCIATED(s1%j%r)) call crap1("EQUAL 4") !call allocw(s1) |
---|
941 | ! call newdacop(S1%j,S2%j) |
---|
942 | ! endif |
---|
943 | END SUBROUTINE EQUAL |
---|
944 | |
---|
945 | SUBROUTINE DEQUAL(R1,S2) |
---|
946 | implicit none |
---|
947 | type (TAYLOR),INTENT(IN)::S2 |
---|
948 | real(dp), INTENT(inOUT)::R1 |
---|
949 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
950 | call check_snake |
---|
951 | |
---|
952 | R1=S2.SUB.'0' |
---|
953 | END SUBROUTINE DEQUAL |
---|
954 | |
---|
955 | SUBROUTINE REQUAL(R1,S2) |
---|
956 | implicit none |
---|
957 | type (TAYLOR),INTENT(IN)::S2 |
---|
958 | REAL(SP), INTENT(inOUT)::R1 |
---|
959 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
960 | |
---|
961 | if(real_warning) call real_stop |
---|
962 | call check_snake |
---|
963 | |
---|
964 | R1=S2.SUB.'0' |
---|
965 | |
---|
966 | END SUBROUTINE REQUAL |
---|
967 | |
---|
968 | function DAABSEQUAL(S2) |
---|
969 | implicit none |
---|
970 | type (TAYLOR),INTENT(IN)::S2 |
---|
971 | real(dp) DAABSEQUAL |
---|
972 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
973 | |
---|
974 | |
---|
975 | DAABSEQUAL=abs(S2.sub.'0') |
---|
976 | |
---|
977 | END function DAABSEQUAL |
---|
978 | |
---|
979 | |
---|
980 | SUBROUTINE DEQUALDACON(S2,R1) |
---|
981 | implicit none |
---|
982 | type (TAYLOR),INTENT(inOUT)::S2 |
---|
983 | real(dp), INTENT(IN)::R1 |
---|
984 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
985 | |
---|
986 | ! if(old) then |
---|
987 | if(s2%i==0) call crap1("DEQUALDACON 1") !call allocw(s2) |
---|
988 | CALL DACON(S2%I,R1) |
---|
989 | ! else |
---|
990 | ! IF (.NOT. ASSOCIATED(s2%j%r)) call crap1("DEQUALDACON 2") !call allocw(s2) |
---|
991 | ! CALL newDACON(S2%j,R1) |
---|
992 | ! endif |
---|
993 | END SUBROUTINE DEQUALDACON |
---|
994 | |
---|
995 | SUBROUTINE EQUALDACON(S2,R1) |
---|
996 | implicit none |
---|
997 | type (TAYLOR),INTENT(inOUT)::S2 |
---|
998 | REAL(SP), INTENT(IN)::R1 |
---|
999 | real(dp) R2 |
---|
1000 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1001 | if(real_warning) call real_stop |
---|
1002 | call check_snake |
---|
1003 | |
---|
1004 | if(real_warning) call real_stop |
---|
1005 | ! if(old) then |
---|
1006 | if(s2%i==0) call crap1("EQUALDACON 1") !call allocw(s2) |
---|
1007 | ! else |
---|
1008 | ! IF (.NOT. ASSOCIATED(s2%j%r)) call crap1("EQUALDACON 2") !call allocw(s2) |
---|
1009 | ! endif |
---|
1010 | r2=REAL(r1,kind=DP) |
---|
1011 | s2=r2 |
---|
1012 | END SUBROUTINE EQUALDACON |
---|
1013 | |
---|
1014 | SUBROUTINE IEQUALDACON(S2,R1) |
---|
1015 | implicit none |
---|
1016 | type (TAYLOR),INTENT(inOUT)::S2 |
---|
1017 | INTEGER, INTENT(IN)::R1 |
---|
1018 | real(dp) r2 |
---|
1019 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1020 | call check_snake |
---|
1021 | |
---|
1022 | |
---|
1023 | ! if(old) then |
---|
1024 | if(s2%i==0) call crap1("IEQUALDACON 1") !call allocw(s2) |
---|
1025 | ! else |
---|
1026 | ! IF (.NOT. ASSOCIATED(s2%j%r)) call crap1("IEQUALDACON 2") !call allocw(s2) |
---|
1027 | ! endif |
---|
1028 | r2=REAL(r1,kind=DP) |
---|
1029 | s2=r2 |
---|
1030 | END SUBROUTINE IEQUALDACON |
---|
1031 | |
---|
1032 | FUNCTION dexpt( S1 ) |
---|
1033 | implicit none |
---|
1034 | TYPE (taylor) dexpt |
---|
1035 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1036 | integer localmaster |
---|
1037 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1038 | localmaster=master |
---|
1039 | |
---|
1040 | ! call check(s1) |
---|
1041 | call ass(dexpt) |
---|
1042 | |
---|
1043 | ! if(old) then |
---|
1044 | call dafun('EXP ',s1%i,temp) |
---|
1045 | call dacop(temp,dexpt%i) |
---|
1046 | ! else |
---|
1047 | ! call newdafun('EXP ',s1%j,dexpt%j) |
---|
1048 | ! endif |
---|
1049 | |
---|
1050 | master=localmaster |
---|
1051 | |
---|
1052 | END FUNCTION dexpt |
---|
1053 | |
---|
1054 | FUNCTION FULL_ABST( S1 ) |
---|
1055 | implicit none |
---|
1056 | real(dp) FULL_ABST |
---|
1057 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1058 | |
---|
1059 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1060 | ! call check(s1) |
---|
1061 | |
---|
1062 | ! if(old) then |
---|
1063 | CALL DAABS(S1%I,FULL_ABST) |
---|
1064 | ! else |
---|
1065 | ! CALL newDAABS(S1%j,FULL_ABST) |
---|
1066 | ! endif |
---|
1067 | |
---|
1068 | END FUNCTION FULL_ABST |
---|
1069 | |
---|
1070 | |
---|
1071 | |
---|
1072 | |
---|
1073 | FUNCTION dtant( S1 ) |
---|
1074 | implicit none |
---|
1075 | TYPE (taylor) dtant |
---|
1076 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1077 | integer localmaster |
---|
1078 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1079 | localmaster=master |
---|
1080 | |
---|
1081 | ! call check(s1) |
---|
1082 | call ass(dtant) |
---|
1083 | |
---|
1084 | ! if(old) then |
---|
1085 | call dafun('SIN ',s1%i,temp) |
---|
1086 | call dacop(temp,dtant%i) |
---|
1087 | call dafun('COS ',s1%i,temp) |
---|
1088 | call dadiv(dtant%i,temp,dtant%i) |
---|
1089 | ! else |
---|
1090 | ! call newdafun('SIN ',s1%j,templ) |
---|
1091 | ! call newdacop(templ,dtant%j) |
---|
1092 | ! call newdafun('COS ',s1%j,templ) |
---|
1093 | ! call newdadiv(dtant%j,templ,dtant%j) |
---|
1094 | ! endif |
---|
1095 | |
---|
1096 | master=localmaster |
---|
1097 | |
---|
1098 | END FUNCTION dtant |
---|
1099 | |
---|
1100 | FUNCTION datanht( S1 ) |
---|
1101 | implicit none |
---|
1102 | TYPE (taylor) datanht |
---|
1103 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1104 | integer localmaster |
---|
1105 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1106 | localmaster=master |
---|
1107 | |
---|
1108 | ! call check(s1) |
---|
1109 | call ass(datanht) |
---|
1110 | |
---|
1111 | datanht=log((1+s1)/sqrt(1-s1))/2.0_dp |
---|
1112 | |
---|
1113 | master=localmaster |
---|
1114 | |
---|
1115 | END FUNCTION datanht |
---|
1116 | |
---|
1117 | FUNCTION dcost( S1 ) |
---|
1118 | implicit none |
---|
1119 | TYPE (taylor) dcost |
---|
1120 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1121 | integer localmaster |
---|
1122 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1123 | localmaster=master |
---|
1124 | |
---|
1125 | |
---|
1126 | |
---|
1127 | ! call check(s1) |
---|
1128 | call ass(dcost) |
---|
1129 | |
---|
1130 | ! if(old) then |
---|
1131 | call dafun('COS ',s1%i,temp) |
---|
1132 | call dacop(temp,dcost%i) |
---|
1133 | ! else |
---|
1134 | ! call newdafun('COS ',s1%j,dcost%j) |
---|
1135 | ! endif |
---|
1136 | |
---|
1137 | master=localmaster |
---|
1138 | |
---|
1139 | END FUNCTION dcost |
---|
1140 | |
---|
1141 | FUNCTION dsint( S1 ) |
---|
1142 | implicit none |
---|
1143 | TYPE (taylor) dsint |
---|
1144 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1145 | integer localmaster |
---|
1146 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1147 | localmaster=master |
---|
1148 | |
---|
1149 | |
---|
1150 | ! call check(s1) |
---|
1151 | call ass(dsint) |
---|
1152 | ! if(old) then |
---|
1153 | call dafun('SIN ',s1%i,temp) |
---|
1154 | call dacop(temp,dsint%i) |
---|
1155 | ! else |
---|
1156 | ! call newdafun('SIN ',s1%j,dsint%j) |
---|
1157 | ! endif |
---|
1158 | |
---|
1159 | master=localmaster |
---|
1160 | |
---|
1161 | END FUNCTION dsint |
---|
1162 | |
---|
1163 | FUNCTION dsinHt( S1 ) |
---|
1164 | implicit none |
---|
1165 | TYPE (taylor) dsinHt |
---|
1166 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1167 | integer localmaster |
---|
1168 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1169 | localmaster=master |
---|
1170 | |
---|
1171 | |
---|
1172 | ! call check(s1) |
---|
1173 | call ass(dsinHt) |
---|
1174 | ! if(old) then |
---|
1175 | call dafun('SINH',s1%i,temp) |
---|
1176 | call dacop(temp,dsinHt%i) |
---|
1177 | ! else |
---|
1178 | ! call newdafun('SINH',s1%j,dsinHt%j) |
---|
1179 | ! endif |
---|
1180 | master=localmaster |
---|
1181 | |
---|
1182 | END FUNCTION dsinHt |
---|
1183 | |
---|
1184 | FUNCTION DCOSHT( S1 ) |
---|
1185 | implicit none |
---|
1186 | TYPE (taylor) DCOSHT |
---|
1187 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1188 | integer localmaster |
---|
1189 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1190 | localmaster=master |
---|
1191 | |
---|
1192 | |
---|
1193 | ! call check(s1) |
---|
1194 | call ass(DCOSHT) |
---|
1195 | ! if(old) then |
---|
1196 | call dafun('COSH',s1%i,temp) |
---|
1197 | call dacop(temp,DCOSHT%i) |
---|
1198 | ! else |
---|
1199 | ! call newdafun('COSH',s1%j,DCOSHT%j) |
---|
1200 | ! endif |
---|
1201 | |
---|
1202 | master=localmaster |
---|
1203 | |
---|
1204 | END FUNCTION DCOSHT |
---|
1205 | |
---|
1206 | FUNCTION dtanht( S1 ) |
---|
1207 | implicit none |
---|
1208 | TYPE (taylor) dtanht |
---|
1209 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1210 | integer localmaster |
---|
1211 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1212 | localmaster=master |
---|
1213 | |
---|
1214 | |
---|
1215 | ! call check(s1) |
---|
1216 | call ass(dtanht) |
---|
1217 | ! if(old) then |
---|
1218 | |
---|
1219 | dtanht=sinh(s1)/cosh(s1) |
---|
1220 | ! else |
---|
1221 | ! call newdafun('COSH',s1%j,DCOSHT%j) |
---|
1222 | ! endif |
---|
1223 | |
---|
1224 | master=localmaster |
---|
1225 | |
---|
1226 | END FUNCTION dtanht |
---|
1227 | |
---|
1228 | FUNCTION dlogt( S1 ) |
---|
1229 | implicit none |
---|
1230 | TYPE (taylor) dlogt |
---|
1231 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1232 | integer localmaster |
---|
1233 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1234 | localmaster=master |
---|
1235 | |
---|
1236 | |
---|
1237 | ! call check(s1) |
---|
1238 | call ass(dlogt) |
---|
1239 | ! if(old) then |
---|
1240 | call dafun('LOG ',s1%i,temp) |
---|
1241 | call dacop(temp,dlogt%i) |
---|
1242 | ! else |
---|
1243 | ! call newdafun('LOG ',s1%j,dlogt%j) |
---|
1244 | ! endif |
---|
1245 | |
---|
1246 | master=localmaster |
---|
1247 | |
---|
1248 | END FUNCTION dlogt |
---|
1249 | |
---|
1250 | FUNCTION dsqrtt( S1 ) |
---|
1251 | implicit none |
---|
1252 | TYPE (taylor) dsqrtt |
---|
1253 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1254 | integer localmaster |
---|
1255 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1256 | localmaster=master |
---|
1257 | |
---|
1258 | |
---|
1259 | ! call check(s1) |
---|
1260 | call ass(dsqrtt) |
---|
1261 | |
---|
1262 | ! if(old) then |
---|
1263 | call dafun('SQRT',s1%i,temp) |
---|
1264 | call dacop(temp,dsqrtt%i) |
---|
1265 | ! else |
---|
1266 | ! call newdafun('SQRT',s1%j,dsqrtt%j) |
---|
1267 | ! endif |
---|
1268 | |
---|
1269 | master=localmaster |
---|
1270 | |
---|
1271 | END FUNCTION dsqrtt |
---|
1272 | |
---|
1273 | FUNCTION mul( S1, S2 ) |
---|
1274 | implicit none |
---|
1275 | TYPE (taylor) mul |
---|
1276 | TYPE (taylor), INTENT (IN) :: S1, S2 |
---|
1277 | integer localmaster |
---|
1278 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1279 | localmaster=master |
---|
1280 | |
---|
1281 | |
---|
1282 | ! call check(s1) |
---|
1283 | ! call check(s2) |
---|
1284 | call ass(mul) |
---|
1285 | |
---|
1286 | ! if(old) then |
---|
1287 | call damul(s1%i,s2%i,temp) |
---|
1288 | call dacop(temp,mul%i) |
---|
1289 | ! else |
---|
1290 | ! call newdamul(s1%j,s2%j,mul%j) |
---|
1291 | ! endif |
---|
1292 | |
---|
1293 | master=localmaster |
---|
1294 | |
---|
1295 | END FUNCTION mul |
---|
1296 | |
---|
1297 | FUNCTION pbbra( S1, S2 ) |
---|
1298 | implicit none |
---|
1299 | TYPE (taylor) pbbra |
---|
1300 | TYPE (taylor), INTENT (IN) :: S1, S2 |
---|
1301 | integer localmaster |
---|
1302 | integer i |
---|
1303 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1304 | localmaster=master |
---|
1305 | |
---|
1306 | |
---|
1307 | ! call check(s1) |
---|
1308 | ! call check(s2) |
---|
1309 | call ass(pbbra) |
---|
1310 | |
---|
1311 | ! if(old) then |
---|
1312 | pbbra=0.0_dp |
---|
1313 | do i=1,nd |
---|
1314 | pbbra=(s1.d.(2*i-1))*(s2.d.(2*i))-(s2.d.(2*i-1))*(s1.d.(2*i))+pbbra |
---|
1315 | enddo |
---|
1316 | ! call DAPOI(s1%i,s2%i,temp,nd) |
---|
1317 | ! call dacop(temp,pbbra%i) |
---|
1318 | ! else |
---|
1319 | ! call newDAPOI(s1%j,s2%j,templ,nd) |
---|
1320 | ! call newdacop(templ,pbbra%j) |
---|
1321 | ! endif |
---|
1322 | |
---|
1323 | master=localmaster |
---|
1324 | |
---|
1325 | END FUNCTION pbbra |
---|
1326 | |
---|
1327 | FUNCTION GETORDER( S1, S2 ) |
---|
1328 | implicit none |
---|
1329 | TYPE (taylor) GETORDER |
---|
1330 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1331 | INTEGER, INTENT (IN) :: S2 |
---|
1332 | integer localmaster |
---|
1333 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1334 | localmaster=master |
---|
1335 | |
---|
1336 | |
---|
1337 | ! call check(s1) |
---|
1338 | call ass(GETORDER) |
---|
1339 | |
---|
1340 | ! if(old) then |
---|
1341 | CALL TAKE(S1%I,S2,TEMP) |
---|
1342 | call dacop(temp,GETORDER%i) |
---|
1343 | ! else |
---|
1344 | ! CALL NEWTAKE(S1%J,S2,TEMPL) |
---|
1345 | ! call NEWdacop(tempL,GETORDER%J) |
---|
1346 | ! endif |
---|
1347 | master=localmaster |
---|
1348 | |
---|
1349 | END FUNCTION GETORDER |
---|
1350 | |
---|
1351 | |
---|
1352 | |
---|
1353 | |
---|
1354 | FUNCTION CUTORDER( S1, S2 ) |
---|
1355 | implicit none |
---|
1356 | TYPE (taylor) CUTORDER |
---|
1357 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1358 | INTEGER, INTENT (IN) :: S2 |
---|
1359 | integer localmaster |
---|
1360 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1361 | localmaster=master |
---|
1362 | |
---|
1363 | ! call check(s1) |
---|
1364 | call ass(CUTORDER) |
---|
1365 | |
---|
1366 | ! if(old) then |
---|
1367 | call datrunc(S1%I,s2,CUTORDER%i) |
---|
1368 | ! call dacop(S1%I,CUTORDER%i) |
---|
1369 | |
---|
1370 | ! DO I=S2,NO |
---|
1371 | ! CALL TAKE(CUTORDER%I,I,TEMP) |
---|
1372 | ! CALL DASUB(CUTORDER%I,TEMP,CUTORDER%I) |
---|
1373 | ! ENDDO |
---|
1374 | ! else |
---|
1375 | ! call NEWdacop(S1%J,CUTORDER%J) |
---|
1376 | ! DO I=S2,NO |
---|
1377 | ! CALL NEWTAKE(CUTORDER%J,I,TEMPL) |
---|
1378 | ! CALL NEWDASUB(CUTORDER%J,TEMPL,CUTORDER%J) |
---|
1379 | ! ENDDO |
---|
1380 | ! endif |
---|
1381 | master=localmaster |
---|
1382 | |
---|
1383 | END FUNCTION CUTORDER |
---|
1384 | |
---|
1385 | FUNCTION dputchar( S1, S2 ) |
---|
1386 | implicit none |
---|
1387 | TYPE (taylor) dputchar |
---|
1388 | real(dp), INTENT (IN) :: S1 |
---|
1389 | CHARACTER(*) , INTENT (IN) :: S2 |
---|
1390 | CHARACTER (LEN = LNV) resul |
---|
1391 | integer j(lnv),i |
---|
1392 | integer localmaster |
---|
1393 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1394 | localmaster=master |
---|
1395 | |
---|
1396 | |
---|
1397 | call ass(dputchar) |
---|
1398 | |
---|
1399 | |
---|
1400 | resul = trim(ADJUSTL (s2)) |
---|
1401 | |
---|
1402 | do i=1,lnv |
---|
1403 | j(i)=0 |
---|
1404 | enddo |
---|
1405 | |
---|
1406 | !frs get around compiler problem |
---|
1407 | nd2par= len(trim(ADJUSTL (s2))) |
---|
1408 | !frs do i=1,len(trim(ADJUSTL (s2))) |
---|
1409 | do i=1,nd2par |
---|
1410 | CALL CHARINT(RESUL(I:I),J(I)) |
---|
1411 | if(i>nv) then |
---|
1412 | if(j(i)>0) then |
---|
1413 | dputchar=0.0_dp |
---|
1414 | ! call var(dputchar,zero,0) |
---|
1415 | return |
---|
1416 | endif |
---|
1417 | endif |
---|
1418 | enddo |
---|
1419 | |
---|
1420 | |
---|
1421 | |
---|
1422 | dputchar=0.0_dp |
---|
1423 | ! call var(dputchar,zero,0) |
---|
1424 | CALL pok(dputchar,j,s1) |
---|
1425 | master=localmaster |
---|
1426 | |
---|
1427 | END FUNCTION dputchar |
---|
1428 | |
---|
1429 | FUNCTION dputint( S1, S2 ) |
---|
1430 | implicit none |
---|
1431 | TYPE (taylor) dputint |
---|
1432 | real(dp), INTENT (IN) :: S1 |
---|
1433 | integer , INTENT (IN) :: S2(:) |
---|
1434 | integer j(lnv),i |
---|
1435 | integer localmaster |
---|
1436 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1437 | localmaster=master |
---|
1438 | |
---|
1439 | |
---|
1440 | call ass(dputint) |
---|
1441 | |
---|
1442 | |
---|
1443 | |
---|
1444 | do i=1,lnv |
---|
1445 | j(i)=0 |
---|
1446 | enddo |
---|
1447 | |
---|
1448 | !frs get around compiler problem |
---|
1449 | nd2par= size(s2) |
---|
1450 | !frs do i=1,len(trim(ADJUSTL (s2))) |
---|
1451 | do i=1,nd2par |
---|
1452 | j(i)=s2(i) |
---|
1453 | enddo |
---|
1454 | do i=1,nd2par |
---|
1455 | if(i>nv) then |
---|
1456 | if(j(i)>0) then |
---|
1457 | ! call var(dputint,zero,0) |
---|
1458 | dputint=0.0_dp |
---|
1459 | return |
---|
1460 | endif |
---|
1461 | endif |
---|
1462 | enddo |
---|
1463 | |
---|
1464 | |
---|
1465 | dputint=0.0_dp |
---|
1466 | ! call var(dputint,zero,0) |
---|
1467 | CALL pok(dputint,j,s1) |
---|
1468 | master=localmaster |
---|
1469 | |
---|
1470 | END FUNCTION dputint |
---|
1471 | |
---|
1472 | FUNCTION dputint0( S1, S2 ) |
---|
1473 | implicit none |
---|
1474 | TYPE (taylor) dputint0 |
---|
1475 | real(dp), INTENT (IN) :: S1 |
---|
1476 | integer , INTENT (IN) :: S2 |
---|
1477 | integer j(lnv) |
---|
1478 | integer localmaster |
---|
1479 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1480 | localmaster=master |
---|
1481 | |
---|
1482 | |
---|
1483 | call ass(dputint0) |
---|
1484 | |
---|
1485 | j=0 |
---|
1486 | if(s2>nv) then |
---|
1487 | dputint0=S1 |
---|
1488 | return |
---|
1489 | endif |
---|
1490 | |
---|
1491 | |
---|
1492 | dputint0=0.0_dp |
---|
1493 | ! call var(dputint0,zero,s2) |
---|
1494 | |
---|
1495 | j(s2)=1 |
---|
1496 | CALL pok(dputint0,j,s1) |
---|
1497 | |
---|
1498 | master=localmaster |
---|
1499 | |
---|
1500 | END FUNCTION dputint0 |
---|
1501 | |
---|
1502 | |
---|
1503 | FUNCTION GETCHARnd2s( S1, S2 ) |
---|
1504 | implicit none |
---|
1505 | TYPE (taylor) GETCHARnd2s |
---|
1506 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1507 | CHARACTER(*) , INTENT (IN) :: S2 |
---|
1508 | |
---|
1509 | integer localmaster |
---|
1510 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1511 | localmaster=master |
---|
1512 | |
---|
1513 | |
---|
1514 | call ass(GETCHARnd2s) |
---|
1515 | |
---|
1516 | |
---|
1517 | GETCHARnd2s=s1.par.s2 |
---|
1518 | call shiftda(GETCHARnd2s,GETCHARnd2s, len(trim(ADJUSTR (s2) ))) |
---|
1519 | |
---|
1520 | master=localmaster |
---|
1521 | |
---|
1522 | |
---|
1523 | END FUNCTION GETCHARnd2s |
---|
1524 | |
---|
1525 | FUNCTION GETintnd2s( S1, S2 ) |
---|
1526 | implicit none |
---|
1527 | TYPE (taylor) GETintnd2s |
---|
1528 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1529 | integer , INTENT (IN) :: S2(:) |
---|
1530 | |
---|
1531 | integer localmaster |
---|
1532 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1533 | localmaster=master |
---|
1534 | |
---|
1535 | |
---|
1536 | call ass(GETintnd2s) |
---|
1537 | |
---|
1538 | |
---|
1539 | GETintnd2s=s1.par.s2 |
---|
1540 | |
---|
1541 | call shiftda(GETintnd2s,GETintnd2s, size(s2) ) |
---|
1542 | |
---|
1543 | master=localmaster |
---|
1544 | |
---|
1545 | |
---|
1546 | END FUNCTION GETintnd2s |
---|
1547 | |
---|
1548 | FUNCTION GETintk( S1, S2 ) |
---|
1549 | implicit none |
---|
1550 | TYPE (taylor) GETintk |
---|
1551 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1552 | integer , INTENT (IN) :: S2 |
---|
1553 | |
---|
1554 | integer localmaster |
---|
1555 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1556 | localmaster=master |
---|
1557 | |
---|
1558 | |
---|
1559 | call ass(GETintk) |
---|
1560 | |
---|
1561 | |
---|
1562 | |
---|
1563 | call shiftda(s1,GETintk, s2 ) |
---|
1564 | |
---|
1565 | master=localmaster |
---|
1566 | |
---|
1567 | |
---|
1568 | END FUNCTION GETintk |
---|
1569 | |
---|
1570 | |
---|
1571 | |
---|
1572 | FUNCTION GETchar( S1, S2 ) |
---|
1573 | |
---|
1574 | implicit none |
---|
1575 | real(dp) GETchar,r1 |
---|
1576 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1577 | CHARACTER(*) , INTENT (IN) :: S2 |
---|
1578 | CHARACTER (LEN = LNV) resul |
---|
1579 | integer j(lnv),i,c |
---|
1580 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1581 | |
---|
1582 | |
---|
1583 | |
---|
1584 | resul = s2 |
---|
1585 | call context(resul) |
---|
1586 | |
---|
1587 | do i=1,lnv |
---|
1588 | j(i)=0 |
---|
1589 | enddo |
---|
1590 | |
---|
1591 | |
---|
1592 | nd2par= len_trim(resul) |
---|
1593 | |
---|
1594 | |
---|
1595 | |
---|
1596 | |
---|
1597 | |
---|
1598 | do i=1,nd2par |
---|
1599 | CALL CHARINT(RESUL(I:I),J(I)) |
---|
1600 | enddo |
---|
1601 | |
---|
1602 | c=0 |
---|
1603 | do i=c_%nv+1,lnv |
---|
1604 | c=j(i)+c |
---|
1605 | enddo |
---|
1606 | |
---|
1607 | if(c>0) then |
---|
1608 | r1=0.0_dp |
---|
1609 | else |
---|
1610 | CALL dapek(S1%I,j,r1) |
---|
1611 | endif |
---|
1612 | |
---|
1613 | |
---|
1614 | GETchar=r1 |
---|
1615 | |
---|
1616 | END FUNCTION GETchar |
---|
1617 | |
---|
1618 | |
---|
1619 | |
---|
1620 | |
---|
1621 | FUNCTION GETint( S1, S2 ) |
---|
1622 | implicit none |
---|
1623 | real(dp) GETint,r1 |
---|
1624 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1625 | integer , INTENT (IN) :: S2(:) |
---|
1626 | integer j(lnv),i,c |
---|
1627 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1628 | |
---|
1629 | |
---|
1630 | |
---|
1631 | do i=1,lnv |
---|
1632 | j(i)=0 |
---|
1633 | enddo |
---|
1634 | |
---|
1635 | |
---|
1636 | nd2par= size(s2) |
---|
1637 | |
---|
1638 | do i=1,nd2par |
---|
1639 | J(I)=s2(i) |
---|
1640 | enddo |
---|
1641 | |
---|
1642 | c=0 |
---|
1643 | do i=c_%nv+1,lnv |
---|
1644 | c=j(i)+c |
---|
1645 | enddo |
---|
1646 | |
---|
1647 | if(c>0) then |
---|
1648 | r1=0.0_dp |
---|
1649 | else |
---|
1650 | CALL dapek(S1%I,j,r1) |
---|
1651 | endif |
---|
1652 | |
---|
1653 | GETint=r1 |
---|
1654 | |
---|
1655 | END FUNCTION GETint |
---|
1656 | |
---|
1657 | |
---|
1658 | |
---|
1659 | |
---|
1660 | FUNCTION GETdiff( S1, S2 ) |
---|
1661 | implicit none |
---|
1662 | TYPE (taylor) GETdiff |
---|
1663 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1664 | INTEGER, INTENT (IN) :: S2 |
---|
1665 | integer localmaster |
---|
1666 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1667 | localmaster=master |
---|
1668 | |
---|
1669 | |
---|
1670 | ! call check(s1) |
---|
1671 | call ass(GETdiff) |
---|
1672 | |
---|
1673 | ! if(old) then |
---|
1674 | CALL dader(S2,S1%I,TEMP) |
---|
1675 | call dacop(temp,GETdiff%i) |
---|
1676 | ! else |
---|
1677 | ! CALL NEWdader(S2,S1%J,TEMPL) |
---|
1678 | ! call NEWdacop(tempL,GETdiff%J) |
---|
1679 | ! endif |
---|
1680 | master=localmaster |
---|
1681 | |
---|
1682 | END FUNCTION GETdiff |
---|
1683 | |
---|
1684 | FUNCTION GETINTegrate( S1, S2 ) |
---|
1685 | implicit none |
---|
1686 | TYPE (taylor) GETINTegrate |
---|
1687 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1688 | INTEGER, INTENT (IN) :: S2 |
---|
1689 | integer localmaster,n,i |
---|
1690 | type(taylor) t,x |
---|
1691 | real(dp) value |
---|
1692 | integer, allocatable :: jc(:) |
---|
1693 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1694 | localmaster=master |
---|
1695 | |
---|
1696 | allocate(jc(c_%nv)) |
---|
1697 | jc=0 |
---|
1698 | ! call check(s1) |
---|
1699 | call ass(GETINTegrate) |
---|
1700 | call alloc(t,x) |
---|
1701 | t=s1 |
---|
1702 | x=0 |
---|
1703 | call taylor_cycle(t,size=n) |
---|
1704 | |
---|
1705 | do i=1,n |
---|
1706 | call taylor_cycle(t,ii=i,value=value,j=jc) |
---|
1707 | |
---|
1708 | x=((value/(jc(s2)+1)).mono.jc)*(1.0_dp.mono.s2)+x |
---|
1709 | |
---|
1710 | enddo |
---|
1711 | |
---|
1712 | GETINTegrate=x |
---|
1713 | |
---|
1714 | call kill(t,x) |
---|
1715 | deallocate(jc) |
---|
1716 | master=localmaster |
---|
1717 | |
---|
1718 | END FUNCTION GETINTegrate |
---|
1719 | |
---|
1720 | |
---|
1721 | FUNCTION GETdatra( S1, S2 ) |
---|
1722 | implicit none |
---|
1723 | TYPE (taylor) GETdatra |
---|
1724 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1725 | INTEGER, INTENT (IN) :: S2 |
---|
1726 | integer localmaster |
---|
1727 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1728 | localmaster=master |
---|
1729 | |
---|
1730 | |
---|
1731 | ! call check(s1) |
---|
1732 | call ass(GETdatra) |
---|
1733 | |
---|
1734 | ! if(old) then |
---|
1735 | CALL datra(S2,S1%I,TEMP) |
---|
1736 | call dacop(temp,GETdatra%i) |
---|
1737 | ! else |
---|
1738 | ! CALL NEWdatra(S2,S1%J,TEMPL) |
---|
1739 | ! call NEWdacop(tempL,GETdatra%J) |
---|
1740 | ! endif |
---|
1741 | master=localmaster |
---|
1742 | |
---|
1743 | END FUNCTION GETdatra |
---|
1744 | |
---|
1745 | FUNCTION POW( S1, R2 ) |
---|
1746 | implicit none |
---|
1747 | TYPE (taylor) POW |
---|
1748 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1749 | INTEGER, INTENT (IN) :: R2 |
---|
1750 | INTEGER I,R22 |
---|
1751 | integer localmaster |
---|
1752 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1753 | localmaster=master |
---|
1754 | |
---|
1755 | |
---|
1756 | ! call check(s1) |
---|
1757 | call ass(POW) |
---|
1758 | |
---|
1759 | ! if(old) then |
---|
1760 | CALL DACON(TEMP,1.0_dp) |
---|
1761 | |
---|
1762 | R22=IABS(R2) |
---|
1763 | DO I=1,R22 |
---|
1764 | CALL DAMUL(TEMP,S1%I,TEMP) |
---|
1765 | ENDDO |
---|
1766 | IF(R2.LT.0) THEN |
---|
1767 | CALL DADIC(TEMP,1.0_dp,TEMP) |
---|
1768 | ENDIF |
---|
1769 | call dacop(temp,POW%i) |
---|
1770 | ! ELSE |
---|
1771 | ! |
---|
1772 | ! CALL newDACON(TEMPl,one) |
---|
1773 | ! |
---|
1774 | ! R22=IABS(R2) |
---|
1775 | ! DO I=1,R22 |
---|
1776 | ! CALL newDAMUL(TEMPl,S1%j,TEMPl) |
---|
1777 | ! ENDDO |
---|
1778 | ! IF(R2.LT.0) THEN |
---|
1779 | ! CALL newDADIC(TEMPl,one,TEMPl) |
---|
1780 | ! ENDIF |
---|
1781 | ! call newdacop(templ,POW%j) |
---|
1782 | ! endif |
---|
1783 | master=localmaster |
---|
1784 | END FUNCTION POW |
---|
1785 | |
---|
1786 | FUNCTION POWR8( S1, R2 ) |
---|
1787 | implicit none |
---|
1788 | TYPE (taylor) POWR8 |
---|
1789 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1790 | real(dp), INTENT (IN) :: R2 |
---|
1791 | integer localmaster |
---|
1792 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1793 | localmaster=master |
---|
1794 | |
---|
1795 | |
---|
1796 | ! call check(s1) |
---|
1797 | call ass(POWR8) |
---|
1798 | |
---|
1799 | ! if(old) then |
---|
1800 | CALL DAFUN('LOG ',S1%I,TEMP) |
---|
1801 | CALL DACMU(TEMP,R2,TEMP) |
---|
1802 | CALL DAFUN('EXP ',TEMP,TEMP) |
---|
1803 | call dacop(temp,POWR8%i) |
---|
1804 | ! ELSE |
---|
1805 | ! CALL NEWDAFUN('LOG ',S1%J,TEMPL) |
---|
1806 | ! CALL NEWDACMU(TEMPL,R2,TEMPL) |
---|
1807 | ! CALL NEWDAFUN('EXP ',TEMPL,POWR8%J) |
---|
1808 | ! endif |
---|
1809 | master=localmaster |
---|
1810 | END FUNCTION POWR8 |
---|
1811 | |
---|
1812 | FUNCTION POWR( S1, R2 ) |
---|
1813 | implicit none |
---|
1814 | TYPE (taylor) POWR |
---|
1815 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1816 | REAL(SP), INTENT (IN) :: R2 |
---|
1817 | integer localmaster |
---|
1818 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1819 | localmaster=master |
---|
1820 | |
---|
1821 | |
---|
1822 | if(real_warning) call real_stop |
---|
1823 | ! call check(s1) |
---|
1824 | call ass(POWR) |
---|
1825 | |
---|
1826 | ! if(old) then |
---|
1827 | CALL DAFUN('LOG ',S1%I,TEMP) |
---|
1828 | CALL DACMU(TEMP,REAL(R2,kind=DP),TEMP) |
---|
1829 | CALL DAFUN('EXP ',TEMP,TEMP) |
---|
1830 | call dacop(temp,POWR%i) |
---|
1831 | ! ELSE |
---|
1832 | ! CALL NEWDAFUN('LOG ',S1%J,TEMPL) |
---|
1833 | ! CALL NEWDACMU(TEMPL,REAL(R2,kind=DP),TEMPL) |
---|
1834 | ! CALL NEWDAFUN('EXP ',TEMPL,POWR%J) |
---|
1835 | ! endif |
---|
1836 | master=localmaster |
---|
1837 | END FUNCTION POWR |
---|
1838 | |
---|
1839 | |
---|
1840 | |
---|
1841 | FUNCTION dmulsc( S1, sc ) |
---|
1842 | implicit none |
---|
1843 | TYPE (taylor) dmulsc |
---|
1844 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1845 | real(dp), INTENT (IN) :: sc |
---|
1846 | integer localmaster |
---|
1847 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1848 | localmaster=master |
---|
1849 | |
---|
1850 | |
---|
1851 | ! call check(s1) |
---|
1852 | call ass(dmulsc) |
---|
1853 | |
---|
1854 | ! if(old) then |
---|
1855 | call dacmu(s1%i,sc,temp) |
---|
1856 | call dacop(temp,dmulsc%i) |
---|
1857 | ! else |
---|
1858 | ! call newdacmu(s1%j,sc,dmulsc%j) |
---|
1859 | ! endif |
---|
1860 | |
---|
1861 | master=localmaster |
---|
1862 | END FUNCTION dmulsc |
---|
1863 | |
---|
1864 | FUNCTION mulsc( S1, sc ) |
---|
1865 | implicit none |
---|
1866 | TYPE (taylor) mulsc |
---|
1867 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1868 | real(sp), INTENT (IN) :: sc |
---|
1869 | integer localmaster |
---|
1870 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1871 | localmaster=master |
---|
1872 | |
---|
1873 | |
---|
1874 | if(real_warning) call real_stop |
---|
1875 | ! call check(s1) |
---|
1876 | call ass(mulsc) |
---|
1877 | |
---|
1878 | ! if(old) then |
---|
1879 | call dacmu(s1%i,REAL(sc,kind=DP),temp) |
---|
1880 | call dacop(temp,mulsc%i) |
---|
1881 | ! else |
---|
1882 | ! call newdacmu(s1%j,REAL(sc,kind=DP),mulsc%j) |
---|
1883 | ! endif |
---|
1884 | master=localmaster |
---|
1885 | END FUNCTION mulsc |
---|
1886 | |
---|
1887 | FUNCTION imulsc( S1, sc ) |
---|
1888 | implicit none |
---|
1889 | TYPE (taylor) imulsc |
---|
1890 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1891 | integer, INTENT (IN) :: sc |
---|
1892 | integer localmaster |
---|
1893 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1894 | localmaster=master |
---|
1895 | |
---|
1896 | |
---|
1897 | ! call check(s1) |
---|
1898 | call ass(imulsc) |
---|
1899 | |
---|
1900 | |
---|
1901 | ! if(old) then |
---|
1902 | call dacmu(s1%i,REAL(sc,kind=DP),temp) |
---|
1903 | call dacop(temp,imulsc%i) |
---|
1904 | ! else |
---|
1905 | ! call newdacmu(s1%j,REAL(sc,kind=DP),imulsc%j) |
---|
1906 | ! endif |
---|
1907 | |
---|
1908 | master=localmaster |
---|
1909 | END FUNCTION imulsc |
---|
1910 | |
---|
1911 | FUNCTION dscmul( sc,S1 ) |
---|
1912 | implicit none |
---|
1913 | TYPE (taylor) dscmul |
---|
1914 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1915 | real(dp), INTENT (IN) :: sc |
---|
1916 | integer localmaster |
---|
1917 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1918 | localmaster=master |
---|
1919 | |
---|
1920 | |
---|
1921 | ! call check(s1) |
---|
1922 | call ass(dscmul) |
---|
1923 | |
---|
1924 | ! if(old) then |
---|
1925 | call dacmu(s1%i,sc,temp) |
---|
1926 | call dacop(temp,dscmul%i) |
---|
1927 | ! else |
---|
1928 | ! call newdacmu(s1%j,sc,dscmul%j) |
---|
1929 | ! endif |
---|
1930 | |
---|
1931 | master=localmaster |
---|
1932 | |
---|
1933 | END FUNCTION dscmul |
---|
1934 | |
---|
1935 | FUNCTION scmul( sc,S1 ) |
---|
1936 | implicit none |
---|
1937 | TYPE (taylor) scmul |
---|
1938 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1939 | real(sp), INTENT (IN) :: sc |
---|
1940 | integer localmaster |
---|
1941 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1942 | localmaster=master |
---|
1943 | |
---|
1944 | |
---|
1945 | if(real_warning) call real_stop |
---|
1946 | ! call check(s1) |
---|
1947 | call ass(scmul) |
---|
1948 | |
---|
1949 | |
---|
1950 | ! if(old) then |
---|
1951 | call dacmu(s1%i,REAL(sc,kind=DP),temp) |
---|
1952 | call dacop(temp,scmul%i) |
---|
1953 | ! else |
---|
1954 | ! call newdacmu(s1%j,REAL(sc,kind=DP),scmul%j) |
---|
1955 | ! endif |
---|
1956 | |
---|
1957 | master=localmaster |
---|
1958 | |
---|
1959 | END FUNCTION scmul |
---|
1960 | |
---|
1961 | FUNCTION iscmul( sc,S1 ) |
---|
1962 | implicit none |
---|
1963 | TYPE (taylor) iscmul |
---|
1964 | TYPE (taylor), INTENT (IN) :: S1 |
---|
1965 | integer, INTENT (IN) :: sc |
---|
1966 | integer localmaster |
---|
1967 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1968 | localmaster=master |
---|
1969 | |
---|
1970 | |
---|
1971 | ! call check(s1) |
---|
1972 | call ass(iscmul) |
---|
1973 | |
---|
1974 | ! if(old) then |
---|
1975 | call dacmu(s1%i,REAL(sc,kind=DP),temp) |
---|
1976 | call dacop(temp,iscmul%i) |
---|
1977 | ! else |
---|
1978 | ! call newdacmu(s1%j,REAL(sc,kind=DP),iscmul%j) |
---|
1979 | ! endif |
---|
1980 | |
---|
1981 | master=localmaster |
---|
1982 | |
---|
1983 | END FUNCTION iscmul |
---|
1984 | |
---|
1985 | FUNCTION div( S1, S2 ) |
---|
1986 | implicit none |
---|
1987 | TYPE (taylor) div |
---|
1988 | TYPE (taylor), INTENT (IN) :: S1, S2 |
---|
1989 | integer localmaster |
---|
1990 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
1991 | localmaster=master |
---|
1992 | |
---|
1993 | |
---|
1994 | ! call check(s1) |
---|
1995 | ! call check(s2) |
---|
1996 | call ass(div) |
---|
1997 | |
---|
1998 | ! if(old) then |
---|
1999 | call dadiv(s1%i,s2%i,temp) |
---|
2000 | call dacop(temp,div%i) |
---|
2001 | ! else |
---|
2002 | ! call newdadiv(s1%j,s2%j,templ) |
---|
2003 | ! call newdacop(templ,div%j) |
---|
2004 | ! endif |
---|
2005 | |
---|
2006 | master=localmaster |
---|
2007 | END FUNCTION div |
---|
2008 | |
---|
2009 | FUNCTION dscdiv( sc,S1 ) |
---|
2010 | implicit none |
---|
2011 | TYPE (taylor) dscdiv |
---|
2012 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2013 | real(dp), INTENT (IN) :: sc |
---|
2014 | integer localmaster |
---|
2015 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2016 | localmaster=master |
---|
2017 | |
---|
2018 | |
---|
2019 | ! call check(s1) |
---|
2020 | call ass(dscdiv) |
---|
2021 | |
---|
2022 | ! if(old) then |
---|
2023 | call dadic(s1%i,sc,temp) |
---|
2024 | call dacop(temp,dscdiv%i) |
---|
2025 | ! else |
---|
2026 | ! call newdadic(s1%j,sc,dscdiv%j) |
---|
2027 | ! endif |
---|
2028 | |
---|
2029 | master=localmaster |
---|
2030 | |
---|
2031 | END FUNCTION dscdiv |
---|
2032 | |
---|
2033 | FUNCTION scdiv( sc,S1 ) |
---|
2034 | implicit none |
---|
2035 | TYPE (taylor) scdiv |
---|
2036 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2037 | real(sp), INTENT (IN) :: sc |
---|
2038 | integer localmaster |
---|
2039 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2040 | localmaster=master |
---|
2041 | |
---|
2042 | |
---|
2043 | if(real_warning) call real_stop |
---|
2044 | ! call check(s1) |
---|
2045 | call ass(scdiv) |
---|
2046 | |
---|
2047 | |
---|
2048 | ! if(old) then |
---|
2049 | call dadic(s1%i,REAL(sc,kind=DP),temp) |
---|
2050 | call dacop(temp,scdiv%i) |
---|
2051 | ! else |
---|
2052 | ! call newdadic(s1%j,REAL(sc,kind=DP),scdiv%j) |
---|
2053 | ! endif |
---|
2054 | |
---|
2055 | master=localmaster |
---|
2056 | END FUNCTION scdiv |
---|
2057 | |
---|
2058 | FUNCTION iscdiv( sc,S1 ) |
---|
2059 | implicit none |
---|
2060 | TYPE (taylor) iscdiv |
---|
2061 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2062 | integer, INTENT (IN) :: sc |
---|
2063 | integer localmaster |
---|
2064 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2065 | localmaster=master |
---|
2066 | |
---|
2067 | |
---|
2068 | ! call check(s1) |
---|
2069 | call ass(iscdiv) |
---|
2070 | |
---|
2071 | ! if(old) then |
---|
2072 | call dadic(s1%i,REAL(sc,kind=DP),temp) |
---|
2073 | call dacop(temp,iscdiv%i) |
---|
2074 | ! else |
---|
2075 | ! call newdadic(s1%j,REAL(sc,kind=DP),iscdiv%j) |
---|
2076 | ! endif |
---|
2077 | |
---|
2078 | |
---|
2079 | master=localmaster |
---|
2080 | END FUNCTION iscdiv |
---|
2081 | |
---|
2082 | FUNCTION ddivsc( S1, sc ) |
---|
2083 | implicit none |
---|
2084 | TYPE (taylor) ddivsc |
---|
2085 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2086 | real(dp), INTENT (IN) :: sc |
---|
2087 | integer localmaster |
---|
2088 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2089 | localmaster=master |
---|
2090 | |
---|
2091 | |
---|
2092 | ! call check(s1) |
---|
2093 | call ass(ddivsc) |
---|
2094 | |
---|
2095 | |
---|
2096 | ! if(old) then |
---|
2097 | call dacdi(s1%i,sc,temp) |
---|
2098 | call dacop(temp,ddivsc%i) |
---|
2099 | ! else |
---|
2100 | ! call newdacdi(s1%j,sc,ddivsc%j) |
---|
2101 | ! endif |
---|
2102 | master=localmaster |
---|
2103 | |
---|
2104 | END FUNCTION ddivsc |
---|
2105 | |
---|
2106 | FUNCTION divsc( S1, sc ) |
---|
2107 | implicit none |
---|
2108 | TYPE (taylor) divsc |
---|
2109 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2110 | real(sp), INTENT (IN) :: sc |
---|
2111 | integer localmaster |
---|
2112 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2113 | localmaster=master |
---|
2114 | |
---|
2115 | |
---|
2116 | if(real_warning) call real_stop |
---|
2117 | ! call check(s1) |
---|
2118 | call ass(divsc) |
---|
2119 | |
---|
2120 | ! if(old) then |
---|
2121 | call dacdi(s1%i,REAL(sc,kind=DP),temp) |
---|
2122 | call dacop(temp,divsc%i) |
---|
2123 | ! else |
---|
2124 | ! call newdacdi(s1%j,REAL(sc,kind=DP),divsc%j) |
---|
2125 | ! endif |
---|
2126 | master=localmaster |
---|
2127 | |
---|
2128 | END FUNCTION divsc |
---|
2129 | |
---|
2130 | |
---|
2131 | FUNCTION idivsc( S1, sc ) |
---|
2132 | implicit none |
---|
2133 | TYPE (taylor) idivsc |
---|
2134 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2135 | integer, INTENT (IN) :: sc |
---|
2136 | integer localmaster |
---|
2137 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2138 | localmaster=master |
---|
2139 | |
---|
2140 | |
---|
2141 | ! call check(s1) |
---|
2142 | call ass(idivsc) |
---|
2143 | |
---|
2144 | |
---|
2145 | ! if(old) then |
---|
2146 | call dacdi(s1%i,REAL(sc,kind=DP),temp) |
---|
2147 | call dacop(temp,idivsc%i) |
---|
2148 | ! else |
---|
2149 | ! call newdacdi(s1%j,REAL(sc,kind=DP),idivsc%j) |
---|
2150 | ! endif |
---|
2151 | master=localmaster |
---|
2152 | |
---|
2153 | END FUNCTION idivsc |
---|
2154 | |
---|
2155 | |
---|
2156 | FUNCTION add( S1, S2 ) |
---|
2157 | implicit none |
---|
2158 | TYPE (taylor) add |
---|
2159 | TYPE (taylor), INTENT (IN) :: S1, S2 |
---|
2160 | integer localmaster |
---|
2161 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2162 | localmaster=master |
---|
2163 | |
---|
2164 | |
---|
2165 | ! call check(s1) |
---|
2166 | ! call check(s2) |
---|
2167 | call ass(add) |
---|
2168 | |
---|
2169 | |
---|
2170 | ! if(old) then |
---|
2171 | call daadd(s1%i,s2%i,add%i) |
---|
2172 | ! call dacop(temp,add%i) |
---|
2173 | ! call daadd(s1%i,s2%i,temp) |
---|
2174 | ! call dacop(temp,add%i) |
---|
2175 | ! else |
---|
2176 | ! call newdaadd(s1%j,s2%j,add%j) |
---|
2177 | ! endif |
---|
2178 | |
---|
2179 | master=localmaster |
---|
2180 | |
---|
2181 | END FUNCTION add |
---|
2182 | |
---|
2183 | FUNCTION daddsc( S1, sc ) |
---|
2184 | implicit none |
---|
2185 | TYPE (taylor) daddsc |
---|
2186 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2187 | real(dp), INTENT (IN) :: sc |
---|
2188 | integer localmaster |
---|
2189 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2190 | localmaster=master |
---|
2191 | |
---|
2192 | |
---|
2193 | ! call check(s1) |
---|
2194 | call ass(daddsc) |
---|
2195 | |
---|
2196 | ! if(old) then |
---|
2197 | call dacad(s1%i,sc,temp) |
---|
2198 | call dacop(temp,daddsc%i) |
---|
2199 | ! else |
---|
2200 | ! call newdacad(s1%j,sc,daddsc%j) |
---|
2201 | ! endif |
---|
2202 | master=localmaster |
---|
2203 | |
---|
2204 | END FUNCTION daddsc |
---|
2205 | |
---|
2206 | FUNCTION addsc( S1, sc ) |
---|
2207 | implicit none |
---|
2208 | TYPE (taylor) addsc |
---|
2209 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2210 | real(sp), INTENT (IN) :: sc |
---|
2211 | integer localmaster |
---|
2212 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2213 | localmaster=master |
---|
2214 | |
---|
2215 | |
---|
2216 | if(real_warning) call real_stop |
---|
2217 | ! call check(s1) |
---|
2218 | call ass(addsc) |
---|
2219 | |
---|
2220 | |
---|
2221 | ! if(old) then |
---|
2222 | call dacad(s1%i,REAL(sc,kind=DP),temp) |
---|
2223 | call dacop(temp,addsc%i) |
---|
2224 | ! else |
---|
2225 | ! call newdacad(s1%j,REAL(sc,kind=DP),addsc%j) |
---|
2226 | ! endif |
---|
2227 | master=localmaster |
---|
2228 | |
---|
2229 | END FUNCTION addsc |
---|
2230 | |
---|
2231 | FUNCTION iaddsc( S1, sc ) |
---|
2232 | implicit none |
---|
2233 | TYPE (taylor) iaddsc |
---|
2234 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2235 | integer, INTENT (IN) :: sc |
---|
2236 | integer localmaster |
---|
2237 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2238 | localmaster=master |
---|
2239 | |
---|
2240 | |
---|
2241 | ! call check(s1) |
---|
2242 | call ass(iaddsc) |
---|
2243 | |
---|
2244 | ! if(old) then |
---|
2245 | call dacad(s1%i,REAL(sc,kind=DP),temp) |
---|
2246 | call dacop(temp,iaddsc%i) |
---|
2247 | ! else |
---|
2248 | ! call newdacad(s1%j,REAL(sc,kind=DP),iaddsc%j) |
---|
2249 | ! endif |
---|
2250 | master=localmaster |
---|
2251 | |
---|
2252 | END FUNCTION iaddsc |
---|
2253 | |
---|
2254 | FUNCTION dscadd( sc,S1) |
---|
2255 | implicit none |
---|
2256 | TYPE (taylor) dscadd |
---|
2257 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2258 | real(dp), INTENT (IN) :: sc |
---|
2259 | integer localmaster |
---|
2260 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2261 | localmaster=master |
---|
2262 | |
---|
2263 | |
---|
2264 | ! call check(s1) |
---|
2265 | call ass(dscadd) |
---|
2266 | |
---|
2267 | ! if(old) then |
---|
2268 | call dacad(s1%i,sc,temp) |
---|
2269 | call dacop(temp,dscadd%i) |
---|
2270 | ! else |
---|
2271 | ! call newdacad(s1%j,sc,dscadd%j) |
---|
2272 | ! endif |
---|
2273 | master=localmaster |
---|
2274 | |
---|
2275 | END FUNCTION dscadd |
---|
2276 | |
---|
2277 | FUNCTION scadd( sc,S1) |
---|
2278 | implicit none |
---|
2279 | TYPE (taylor) scadd |
---|
2280 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2281 | real(sp), INTENT (IN) :: sc |
---|
2282 | integer localmaster |
---|
2283 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2284 | localmaster=master |
---|
2285 | |
---|
2286 | |
---|
2287 | if(real_warning) call real_stop |
---|
2288 | ! call check(s1) |
---|
2289 | call ass(scadd) |
---|
2290 | |
---|
2291 | ! if(old) then |
---|
2292 | call dacad(s1%i,REAL(sc,kind=DP),temp) |
---|
2293 | call dacop(temp,scadd%i) |
---|
2294 | ! else |
---|
2295 | ! call newdacad(s1%j,REAL(sc,kind=DP),scadd%j) |
---|
2296 | ! endif |
---|
2297 | |
---|
2298 | master=localmaster |
---|
2299 | |
---|
2300 | END FUNCTION scadd |
---|
2301 | |
---|
2302 | FUNCTION iscadd( sc,S1) |
---|
2303 | implicit none |
---|
2304 | TYPE (taylor) iscadd |
---|
2305 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2306 | integer, INTENT (IN) :: sc |
---|
2307 | integer localmaster |
---|
2308 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2309 | localmaster=master |
---|
2310 | |
---|
2311 | |
---|
2312 | ! call check(s1) |
---|
2313 | call ass(iscadd) |
---|
2314 | |
---|
2315 | |
---|
2316 | ! if(old) then |
---|
2317 | call dacad(s1%i,REAL(sc,kind=DP),temp) |
---|
2318 | call dacop(temp,iscadd%i) |
---|
2319 | ! else |
---|
2320 | ! call newdacad(s1%j,REAL(sc,kind=DP),iscadd%j) |
---|
2321 | ! endif |
---|
2322 | master=localmaster |
---|
2323 | |
---|
2324 | END FUNCTION iscadd |
---|
2325 | |
---|
2326 | FUNCTION subs( S1, S2 ) |
---|
2327 | implicit none |
---|
2328 | TYPE (taylor) subs |
---|
2329 | TYPE (taylor), INTENT (IN) :: S1, S2 |
---|
2330 | integer localmaster |
---|
2331 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2332 | localmaster=master |
---|
2333 | |
---|
2334 | |
---|
2335 | ! call check(s1) |
---|
2336 | ! call check(s2) |
---|
2337 | call ass(subs) |
---|
2338 | |
---|
2339 | |
---|
2340 | |
---|
2341 | ! if(old) then |
---|
2342 | call dasub(s1%i,s2%i,temp) |
---|
2343 | call dacop(temp,subs%i) |
---|
2344 | ! else |
---|
2345 | ! call newdasub(s1%j,s2%j,subs%j) |
---|
2346 | ! endif |
---|
2347 | master=localmaster |
---|
2348 | |
---|
2349 | END FUNCTION subs |
---|
2350 | |
---|
2351 | FUNCTION dsubsc( S1, sc ) |
---|
2352 | implicit none |
---|
2353 | TYPE (taylor) dsubsc |
---|
2354 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2355 | real(dp), INTENT (IN) :: sc |
---|
2356 | integer localmaster |
---|
2357 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2358 | localmaster=master |
---|
2359 | |
---|
2360 | |
---|
2361 | ! call check(s1) |
---|
2362 | call ass(dsubsc) |
---|
2363 | |
---|
2364 | ! if(old) then |
---|
2365 | call dacsu(s1%i,sc,temp) |
---|
2366 | call dacop(temp,dsubsc%i) |
---|
2367 | ! else |
---|
2368 | ! call newdacsu(s1%j,sc,dsubsc%j) |
---|
2369 | ! endif |
---|
2370 | |
---|
2371 | master=localmaster |
---|
2372 | |
---|
2373 | |
---|
2374 | END FUNCTION dsubsc |
---|
2375 | |
---|
2376 | FUNCTION subsc( S1, sc ) |
---|
2377 | implicit none |
---|
2378 | TYPE (taylor) subsc |
---|
2379 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2380 | real(sp), INTENT (IN) :: sc |
---|
2381 | integer localmaster |
---|
2382 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2383 | localmaster=master |
---|
2384 | |
---|
2385 | |
---|
2386 | if(real_warning) call real_stop |
---|
2387 | ! call check(s1) |
---|
2388 | call ass(subsc) |
---|
2389 | |
---|
2390 | ! if(old) then |
---|
2391 | call dacsu(s1%i,REAL(sc,kind=DP),temp) |
---|
2392 | call dacop(temp,subsc%i) |
---|
2393 | ! else |
---|
2394 | ! call newdacsu(s1%j,REAL(sc,kind=DP),subsc%j) |
---|
2395 | ! endif |
---|
2396 | master=localmaster |
---|
2397 | |
---|
2398 | END FUNCTION subsc |
---|
2399 | |
---|
2400 | FUNCTION isubsc( S1, sc ) |
---|
2401 | implicit none |
---|
2402 | TYPE (taylor) isubsc |
---|
2403 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2404 | integer, INTENT (IN) :: sc |
---|
2405 | integer localmaster |
---|
2406 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2407 | localmaster=master |
---|
2408 | |
---|
2409 | |
---|
2410 | ! call check(s1) |
---|
2411 | call ass(isubsc) |
---|
2412 | |
---|
2413 | ! if(old) then |
---|
2414 | call dacsu(s1%i,REAL(sc,kind=DP),temp) |
---|
2415 | call dacop(temp,isubsc%i) |
---|
2416 | ! else |
---|
2417 | ! call newdacsu(s1%j,REAL(sc,kind=DP),isubsc%j) |
---|
2418 | ! endif |
---|
2419 | master=localmaster |
---|
2420 | |
---|
2421 | END FUNCTION isubsc |
---|
2422 | |
---|
2423 | FUNCTION dscsub( sc,S1) |
---|
2424 | implicit none |
---|
2425 | TYPE (taylor) dscsub |
---|
2426 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2427 | real(dp), INTENT (IN) :: sc |
---|
2428 | integer localmaster |
---|
2429 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2430 | localmaster=master |
---|
2431 | |
---|
2432 | |
---|
2433 | ! call check(s1) |
---|
2434 | call ass(dscsub) |
---|
2435 | |
---|
2436 | ! if(old) then |
---|
2437 | call dasuc(s1%i,sc,temp) |
---|
2438 | call dacop(temp,dscsub%i) |
---|
2439 | ! else |
---|
2440 | ! call newdasuc(s1%j,sc,dscsub%j) |
---|
2441 | ! endif |
---|
2442 | master=localmaster |
---|
2443 | |
---|
2444 | END FUNCTION dscsub |
---|
2445 | |
---|
2446 | FUNCTION scsub( sc,S1) |
---|
2447 | implicit none |
---|
2448 | TYPE (taylor) scsub |
---|
2449 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2450 | real(sp), INTENT (IN) :: sc |
---|
2451 | integer localmaster |
---|
2452 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2453 | localmaster=master |
---|
2454 | |
---|
2455 | |
---|
2456 | if(real_warning) call real_stop |
---|
2457 | ! call check(s1) |
---|
2458 | call ass(scsub) |
---|
2459 | |
---|
2460 | ! if(old) then |
---|
2461 | call dasuc(s1%i,REAL(sc,kind=DP),temp) |
---|
2462 | call dacop(temp,scsub%i) |
---|
2463 | ! else |
---|
2464 | ! call newdasuc(s1%j,REAL(sc,kind=DP),scsub%j) |
---|
2465 | ! endif |
---|
2466 | master=localmaster |
---|
2467 | |
---|
2468 | END FUNCTION scsub |
---|
2469 | |
---|
2470 | FUNCTION iscsub( sc,S1) |
---|
2471 | implicit none |
---|
2472 | TYPE (taylor) iscsub |
---|
2473 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2474 | integer, INTENT (IN) :: sc |
---|
2475 | integer localmaster |
---|
2476 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2477 | localmaster=master |
---|
2478 | |
---|
2479 | |
---|
2480 | ! call check(s1) |
---|
2481 | call ass(iscsub) |
---|
2482 | |
---|
2483 | ! if(old) then |
---|
2484 | call dasuc(s1%i,REAL(sc,kind=DP),temp) |
---|
2485 | call dacop(temp,iscsub%i) |
---|
2486 | ! else |
---|
2487 | ! call newdasuc(s1%j,REAL(sc,kind=DP),iscsub%j) |
---|
2488 | ! endif |
---|
2489 | master=localmaster |
---|
2490 | |
---|
2491 | END FUNCTION iscsub |
---|
2492 | |
---|
2493 | |
---|
2494 | ! These are new general TPSA-Routines |
---|
2495 | |
---|
2496 | |
---|
2497 | |
---|
2498 | FUNCTION varf( S1, S2 ) |
---|
2499 | implicit none |
---|
2500 | TYPE (taylor) varf |
---|
2501 | real(dp), INTENT (IN) :: S1 |
---|
2502 | integer , INTENT (IN) :: S2 |
---|
2503 | integer localmaster |
---|
2504 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2505 | localmaster=master |
---|
2506 | |
---|
2507 | |
---|
2508 | call ass(varf) |
---|
2509 | |
---|
2510 | varf=S1 + (1.0_dp.mono.S2) |
---|
2511 | |
---|
2512 | master=localmaster |
---|
2513 | |
---|
2514 | END FUNCTION varf |
---|
2515 | |
---|
2516 | FUNCTION varf001( S1, S2 ) |
---|
2517 | implicit none |
---|
2518 | TYPE (taylor) varf001 |
---|
2519 | real(dp), INTENT (IN) :: S1(2) |
---|
2520 | integer , INTENT (IN) :: S2 |
---|
2521 | integer localmaster |
---|
2522 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2523 | localmaster=master |
---|
2524 | |
---|
2525 | |
---|
2526 | call ass(varf001) |
---|
2527 | |
---|
2528 | varf001=S1(1) + (s1(2).mono.S2) |
---|
2529 | |
---|
2530 | master=localmaster |
---|
2531 | |
---|
2532 | END FUNCTION varf001 |
---|
2533 | |
---|
2534 | |
---|
2535 | |
---|
2536 | SUBROUTINE shift000(S1,S2,s) |
---|
2537 | implicit none |
---|
2538 | INTEGER,INTENT(IN)::s |
---|
2539 | type (taylor),INTENT(IN)::S1 |
---|
2540 | type (taylor),INTENT(inout)::S2 |
---|
2541 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2542 | |
---|
2543 | ! if(old) then |
---|
2544 | if(s2%i==0) call crap1("shift000 1" ) !call etall1(s2%i) |
---|
2545 | CALL DAshift(s1%i,s2%i,s) |
---|
2546 | ! else |
---|
2547 | ! if(.NOT. ASSOCIATED(s2%j%r))call crap1("shift000 2" ) ! call newetall(s2%j,1) |
---|
2548 | ! |
---|
2549 | ! CALL NEWDAshift(s1%j,s2%j,s) |
---|
2550 | ! endif |
---|
2551 | ! |
---|
2552 | END SUBROUTINE shift000 |
---|
2553 | |
---|
2554 | |
---|
2555 | SUBROUTINE pek000(S1,J,R1) |
---|
2556 | implicit none |
---|
2557 | INTEGER,INTENT(IN),dimension(:)::j |
---|
2558 | real(dp),INTENT(inOUT)::R1 |
---|
2559 | type (taylor),INTENT(IN)::S1 |
---|
2560 | ! integer k |
---|
2561 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2562 | ! if(old) then |
---|
2563 | if(s1%i==0) call crap1("pek000 1" ) !call etall1(s1%i) |
---|
2564 | ! k=s1%i |
---|
2565 | ! write(6,*) r1,k |
---|
2566 | CALL DApek(s1%i,j,r1) |
---|
2567 | ! else |
---|
2568 | ! if(.NOT. ASSOCIATED(s1%j%r)) call crap1("pek000 2" ) ! newetall(s1%j,1) |
---|
2569 | ! |
---|
2570 | ! CALL newDApek(s1%j,j,r1) |
---|
2571 | ! endif |
---|
2572 | ! |
---|
2573 | END SUBROUTINE pek000 |
---|
2574 | |
---|
2575 | SUBROUTINE pok000(S1,J,R1) |
---|
2576 | implicit none |
---|
2577 | INTEGER,INTENT(in),dimension(:)::j |
---|
2578 | real(dp),INTENT(in)::R1 |
---|
2579 | type (taylor),INTENT(inout)::S1 |
---|
2580 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2581 | |
---|
2582 | if(check_j(j)/=0) return |
---|
2583 | ! if(old) then |
---|
2584 | if(s1%i==0) call crap1("pok000 1" ) ! call etall1(s1%i) |
---|
2585 | CALL DApok(s1%i,j,r1) |
---|
2586 | ! else |
---|
2587 | ! if(.NOT. ASSOCIATED(s1%j%r)) call crap1("pok000 2" ) ! call newetall(s1%j,1) |
---|
2588 | ! |
---|
2589 | ! CALL newDApok(s1%j,j,r1) |
---|
2590 | ! endif |
---|
2591 | ! |
---|
2592 | END SUBROUTINE pok000 |
---|
2593 | |
---|
2594 | SUBROUTINE TAYLOR_ran(S1,r1,R2) |
---|
2595 | implicit none |
---|
2596 | real(dp),INTENT(in)::R1 |
---|
2597 | real(dp),INTENT(inout)::R2 |
---|
2598 | type (taylor),INTENT(inout)::S1 |
---|
2599 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2600 | |
---|
2601 | ! |
---|
2602 | ! THIS SUBROUTINE FILLS THE DA VECTOR A WITH RANDOM ENTRIES. |
---|
2603 | ! FOR R1 > 0, THE VECTOR IS FILLED WITH REALS, |
---|
2604 | ! FOR R1 < 0, THE VECTOR IS FILLED WITH SINGLE DIGIT INTEGERS |
---|
2605 | ! ABS(R1) IS THE FILLING FACTOR |
---|
2606 | ! if(old) then |
---|
2607 | if(s1%i==0) call crap1("tAYLOR_ran 1" ) ! call etall1(s1%i) |
---|
2608 | call daran(s1%i,r1,R2) |
---|
2609 | ! else |
---|
2610 | ! if(.NOT. ASSOCIATED(s1%j%r))call crap1("tAYLOR_ran 2" ) ! call newetall(s1%j,1) |
---|
2611 | ! |
---|
2612 | ! call newdaran(s1%j,r1,R2) |
---|
2613 | ! endif |
---|
2614 | ! |
---|
2615 | END SUBROUTINE TAYLOR_ran |
---|
2616 | |
---|
2617 | SUBROUTINE intd_taylor(S1,S2,factor) |
---|
2618 | implicit none |
---|
2619 | type (taylor),INTENT(inOUT)::S2 |
---|
2620 | type (taylor),INTENT(IN)::S1(:) |
---|
2621 | real(dp),INTENT(IN):: factor |
---|
2622 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2623 | |
---|
2624 | ! if(old) then |
---|
2625 | if(s1(1)%i==0) call crap1("intd_taylor 1") !call etall1(s2%h%i) |
---|
2626 | CALL intd(S1%i,s2%i,factor) |
---|
2627 | ! else |
---|
2628 | ! if(.NOT. ASSOCIATED(s1(1)%j%r)) call crap1("intd_taylor 2") !call etall1(s2%h%i) |
---|
2629 | ! CALL newintd(S1%j,s2%j,factor) |
---|
2630 | ! endif |
---|
2631 | END SUBROUTINE intd_taylor |
---|
2632 | |
---|
2633 | SUBROUTINE DIFd_taylor(S2,S1,factor) |
---|
2634 | implicit none |
---|
2635 | type (taylor),INTENT(in)::S2 |
---|
2636 | type (taylor),INTENT(INOUT)::S1(:) |
---|
2637 | real(dp),INTENT(IN):: factor |
---|
2638 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2639 | ! if(old) then |
---|
2640 | CALL DIFD(S2%i,s1%i,factor) |
---|
2641 | ! else |
---|
2642 | ! CALL NEWDIFD(S2%j,s1%j,factor) |
---|
2643 | ! endif |
---|
2644 | END SUBROUTINE DIFd_taylor |
---|
2645 | |
---|
2646 | |
---|
2647 | SUBROUTINE CFU000(S2,FUN,S1) |
---|
2648 | implicit none |
---|
2649 | type (taylor),INTENT(INOUT)::S1 |
---|
2650 | type (taylor),INTENT(IN)::S2 |
---|
2651 | real(dp) FUN |
---|
2652 | EXTERNAL FUN |
---|
2653 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2654 | ! if(old) then |
---|
2655 | if(s1%i==0) call crap1("CFU000 1" ) ! call etall1(s1%i) |
---|
2656 | CALL DACFU(s2%i,FUN,s1%i) |
---|
2657 | ! else |
---|
2658 | ! if(.NOT. ASSOCIATED(s1%j%r))call crap1("CFU000 2" ) ! call newetall(s1%j,1) |
---|
2659 | ! CALL NEWDACFU(s2%J,FUN,s1%J) |
---|
2660 | ! endif |
---|
2661 | |
---|
2662 | END SUBROUTINE CFU000 |
---|
2663 | |
---|
2664 | SUBROUTINE CFUR(S2,FUN,S1) |
---|
2665 | implicit none |
---|
2666 | type (taylor),INTENT(INOUT)::S1 |
---|
2667 | type (taylor),INTENT(IN)::S2 |
---|
2668 | complex(dp) FUN |
---|
2669 | EXTERNAL FUN |
---|
2670 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2671 | ! if(old) then |
---|
2672 | if(s1%i==0) call crap1("CFUR 1" ) ! call etall1(s1%i) |
---|
2673 | CALL DACFUR(s2%i,FUN,s1%i) |
---|
2674 | ! else |
---|
2675 | ! if(.NOT. ASSOCIATED(s1%j%r))call crap1("CFUR 2" ) ! call newetall(s1%j,1) |
---|
2676 | ! CALL NEWDACFUR(s2%J,FUN,s1%J) |
---|
2677 | ! endif |
---|
2678 | |
---|
2679 | END SUBROUTINE CFUR |
---|
2680 | |
---|
2681 | SUBROUTINE CFUI(S2,FUN,S1) |
---|
2682 | implicit none |
---|
2683 | type (taylor),INTENT(INOUT)::S1 |
---|
2684 | type (taylor),INTENT(IN)::S2 |
---|
2685 | complex(dp) FUN |
---|
2686 | EXTERNAL FUN |
---|
2687 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2688 | ! if(old) then |
---|
2689 | if(s1%i==0)call crap1("CFUI 1" ) ! call etall1(s1%i) |
---|
2690 | CALL DACFUI(s2%i,FUN,s1%i) |
---|
2691 | ! else |
---|
2692 | ! if(.NOT. ASSOCIATED(s1%j%r)) call crap1("CFUI 2" ) !call newetall(s1%j,1) |
---|
2693 | ! CALL NEWDACFUI(s2%J,FUN,s1%J) |
---|
2694 | ! endif |
---|
2695 | |
---|
2696 | END SUBROUTINE CFUI |
---|
2697 | |
---|
2698 | SUBROUTINE taylor_eps(r1) |
---|
2699 | implicit none |
---|
2700 | real(dp),INTENT(INOUT)::r1 |
---|
2701 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2702 | ! if(old) then |
---|
2703 | CALL DAeps(r1) |
---|
2704 | ! else |
---|
2705 | ! CALL newDAeps(r1) |
---|
2706 | ! endif |
---|
2707 | |
---|
2708 | END SUBROUTINE taylor_eps |
---|
2709 | |
---|
2710 | |
---|
2711 | |
---|
2712 | FUNCTION GETCHARnd2( S1, S2 ) |
---|
2713 | implicit none |
---|
2714 | TYPE (taylor) GETCHARnd2,junk |
---|
2715 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2716 | CHARACTER(*) , INTENT (IN) :: S2 |
---|
2717 | CHARACTER (LEN = LNV) resul |
---|
2718 | integer i,k |
---|
2719 | integer localmaster |
---|
2720 | |
---|
2721 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2722 | localmaster=master |
---|
2723 | |
---|
2724 | ndel=0 |
---|
2725 | ! call check(s1) |
---|
2726 | call ass(GETCHARnd2) |
---|
2727 | |
---|
2728 | call alloc(junk) |
---|
2729 | resul = trim(ADJUSTR (s2)) |
---|
2730 | |
---|
2731 | do i=1,lnv |
---|
2732 | jfil(i)=0 |
---|
2733 | enddo |
---|
2734 | |
---|
2735 | nd2par= len(trim(ADJUSTR (s2))) |
---|
2736 | |
---|
2737 | !frs get around compiler problem |
---|
2738 | !frs do i=1,len(trim(ADJUSTR (s2))) |
---|
2739 | do i=1,nd2par |
---|
2740 | CALL CHARINT(RESUL(I:I),Jfil(I)) |
---|
2741 | if(i>nv) then |
---|
2742 | if(Jfil(i)>0) then |
---|
2743 | GETCHARnd2=0.0_dp |
---|
2744 | return |
---|
2745 | endif |
---|
2746 | endif |
---|
2747 | |
---|
2748 | enddo |
---|
2749 | |
---|
2750 | !do i=nd2+ndel+1,nv |
---|
2751 | do i=nd2par+1,nv |
---|
2752 | if(jfil(i)/=0) then |
---|
2753 | w_p=0 |
---|
2754 | w_p%nc=1 |
---|
2755 | w_p%fc='(1((1X,A72),/))' |
---|
2756 | w_p%c(1)=" error in getchar for .para. " |
---|
2757 | ! call !write_e(0) |
---|
2758 | stop |
---|
2759 | endif |
---|
2760 | enddo |
---|
2761 | |
---|
2762 | call cfu(s1,filter,junk) |
---|
2763 | |
---|
2764 | !DO I=1,ND2+ndel |
---|
2765 | DO I=1,ND2par |
---|
2766 | DO K=1,JFIL(I) |
---|
2767 | JUNK=JUNK.K.I |
---|
2768 | ENDDO |
---|
2769 | ENDDO |
---|
2770 | |
---|
2771 | GETCHARnd2=junk |
---|
2772 | |
---|
2773 | call kill(junk) |
---|
2774 | master=localmaster |
---|
2775 | |
---|
2776 | END FUNCTION GETCHARnd2 |
---|
2777 | |
---|
2778 | FUNCTION GETintnd2( S1, S2 ) |
---|
2779 | implicit none |
---|
2780 | TYPE (taylor) GETintnd2,junk |
---|
2781 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2782 | integer , INTENT (IN) :: S2(:) |
---|
2783 | integer i,k |
---|
2784 | integer localmaster |
---|
2785 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2786 | localmaster=master |
---|
2787 | |
---|
2788 | ! call check(s1) |
---|
2789 | call ass(GETintnd2) |
---|
2790 | |
---|
2791 | call alloc(junk) |
---|
2792 | |
---|
2793 | do i=1,lnv |
---|
2794 | jfil(i)=0 |
---|
2795 | enddo |
---|
2796 | nd2par=size(s2) |
---|
2797 | ndel=0 |
---|
2798 | |
---|
2799 | !frs get around compiler problem |
---|
2800 | !frs do i=1,len(trim(ADJUSTR (s2))) |
---|
2801 | do i=1,nd2par |
---|
2802 | Jfil(I)=s2(i) |
---|
2803 | if(i>nv) then |
---|
2804 | if(Jfil(i)>0) then |
---|
2805 | GETintnd2=0.0_dp |
---|
2806 | return |
---|
2807 | endif |
---|
2808 | endif |
---|
2809 | |
---|
2810 | enddo |
---|
2811 | |
---|
2812 | !do i=nd2+ndel+1,nv |
---|
2813 | do i=nd2par+1,nv |
---|
2814 | if(jfil(i)/=0) then |
---|
2815 | w_p=0 |
---|
2816 | w_p%nc=1 |
---|
2817 | w_p%fc='(1((1X,A72),/))' |
---|
2818 | w_p%c(1)=" error in GETintnd2 for .para. " |
---|
2819 | ! call !write_e(0) |
---|
2820 | stop |
---|
2821 | endif |
---|
2822 | enddo |
---|
2823 | |
---|
2824 | call cfu(s1,filter,junk) |
---|
2825 | |
---|
2826 | !DO I=1,ND2+ndel |
---|
2827 | DO I=1,ND2par |
---|
2828 | DO K=1,JFIL(I) |
---|
2829 | JUNK=JUNK.K.I |
---|
2830 | ENDDO |
---|
2831 | ENDDO |
---|
2832 | |
---|
2833 | GETintnd2=junk |
---|
2834 | |
---|
2835 | call kill(junk) |
---|
2836 | master=localmaster |
---|
2837 | |
---|
2838 | END FUNCTION GETintnd2 |
---|
2839 | |
---|
2840 | FUNCTION GETintnd2t( S1, S22 ) |
---|
2841 | implicit none |
---|
2842 | TYPE (taylor) GETintnd2t,junk |
---|
2843 | TYPE (taylor), INTENT (IN) :: S1 |
---|
2844 | type(sub_taylor), INTENT (IN) :: S22 |
---|
2845 | integer s2(lnv) |
---|
2846 | integer i |
---|
2847 | integer localmaster |
---|
2848 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2849 | localmaster=master |
---|
2850 | |
---|
2851 | ! call check(s1) |
---|
2852 | call ass(GETintnd2t) |
---|
2853 | |
---|
2854 | call alloc(junk) |
---|
2855 | |
---|
2856 | do i=1,lnv |
---|
2857 | jfilt(i)=0 |
---|
2858 | enddo |
---|
2859 | s2=s22%j |
---|
2860 | nd2part=s22%min |
---|
2861 | nd2partt=s22%max |
---|
2862 | ndel=0 |
---|
2863 | !frs get around compiler problem |
---|
2864 | !frs do i=1,len(trim(ADJUSTR (s2))) |
---|
2865 | do i=nd2part,nd2partt |
---|
2866 | jfilt(I)=s2(i) |
---|
2867 | if(i>nv) then |
---|
2868 | if(jfilt(i)>0) then |
---|
2869 | GETintnd2t=0.0_dp |
---|
2870 | return |
---|
2871 | endif |
---|
2872 | endif |
---|
2873 | |
---|
2874 | enddo |
---|
2875 | |
---|
2876 | !do i=nd2+ndel+1,nv |
---|
2877 | do i=nd2partt+1,nv |
---|
2878 | if(jfilt(i)/=0) then |
---|
2879 | w_p=0 |
---|
2880 | w_p%nc=1 |
---|
2881 | w_p%fc='(1((1X,A72),/))' |
---|
2882 | w_p%c(1)=" error in GETintnd2t for .part_taylor. " |
---|
2883 | ! call !write_e(0) |
---|
2884 | stop |
---|
2885 | endif |
---|
2886 | enddo |
---|
2887 | |
---|
2888 | call cfu(s1,filter_part,junk) |
---|
2889 | |
---|
2890 | !DO I=1,ND2+ndel |
---|
2891 | ! DO I=1,ND2par |
---|
2892 | ! DO K=1,jfilt(I) |
---|
2893 | ! JUNK=JUNK.K.I |
---|
2894 | ! ENDDO |
---|
2895 | ! ENDDO |
---|
2896 | |
---|
2897 | GETintnd2t=junk |
---|
2898 | |
---|
2899 | call kill(junk) |
---|
2900 | master=localmaster |
---|
2901 | |
---|
2902 | END FUNCTION GETintnd2t |
---|
2903 | |
---|
2904 | |
---|
2905 | SUBROUTINE taylor_cycle(S1,size,ii,VALUE,J) |
---|
2906 | implicit none |
---|
2907 | type (taylor),INTENT(IN)::S1 |
---|
2908 | integer,optional, intent(inout):: size |
---|
2909 | integer,optional, intent(in):: ii |
---|
2910 | integer,optional, intent(inout)::J(:) |
---|
2911 | real(dp), OPTIONAL, intent(inout):: value |
---|
2912 | INTEGER ipresent,ILLA |
---|
2913 | real(dp) VALUE0 |
---|
2914 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
2915 | ! if(old) THEN |
---|
2916 | IF(PRESENT(J).AND.PRESENT(VALUE).and.present(ii)) THEN |
---|
2917 | call dacycle(S1%i,ii,value,illa,J) |
---|
2918 | ELSEif(present(size)) then |
---|
2919 | call dacycle(S1%i,ipresent,value0,size) |
---|
2920 | else |
---|
2921 | write(6,*) "error in taylor_cycle" |
---|
2922 | stop 888 |
---|
2923 | ENDIF |
---|
2924 | |
---|
2925 | END SUBROUTINE taylor_cycle |
---|
2926 | |
---|
2927 | |
---|
2928 | ! SUBROUTINE taylor_clean(S1,VALUE) |
---|
2929 | ! implicit none |
---|
2930 | ! type (taylor),INTENT(INout)::S1 |
---|
2931 | ! real(dp) value |
---|
2932 | ! call daclean(S1%i,value) |
---|
2933 | ! END SUBROUTINE taylor_clean |
---|
2934 | |
---|
2935 | subroutine check_snake() |
---|
2936 | implicit none |
---|
2937 | master=master+1 |
---|
2938 | select case (master) |
---|
2939 | case(1:ndumt) |
---|
2940 | if(iass0user(master)>scratchda(master)%n.or.scratchda(master)%n>newscheme_max) then |
---|
2941 | w_p=0 |
---|
2942 | w_p%nc=1 |
---|
2943 | w_p%fc='(1((1X,A72),/))' |
---|
2944 | w_p%fi='(3((1X,i4)))' |
---|
2945 | w_p%c(1)= "iass0user(master),scratchda(master)%n,newscheme_max" |
---|
2946 | w_p=(/iass0user(master),scratchda(master)%n,newscheme_max/) |
---|
2947 | ! call !write_e |
---|
2948 | call ndum_warning_user |
---|
2949 | endif |
---|
2950 | iass0user(master)=0 |
---|
2951 | case(ndumt+1:) |
---|
2952 | w_p=0 |
---|
2953 | w_p%nc=1 |
---|
2954 | w_p=(/"Should not be here"/) |
---|
2955 | w_p%fc='(1((1X,A72),/))' |
---|
2956 | ! call !write_e(101) |
---|
2957 | end select |
---|
2958 | master=master-1 |
---|
2959 | end subroutine check_snake |
---|
2960 | |
---|
2961 | ! functions used inside other routines |
---|
2962 | |
---|
2963 | SUBROUTINE CHARINT(A,I) |
---|
2964 | IMPLICIT NONE |
---|
2965 | INTEGER I |
---|
2966 | CHARACTER(1) A |
---|
2967 | |
---|
2968 | i=-1 |
---|
2969 | IF(A=='1') I=1 |
---|
2970 | IF(A=='2') I=2 |
---|
2971 | IF(A=='3') I=3 |
---|
2972 | IF(A=='4') I=4 |
---|
2973 | IF(A=='5') I=5 |
---|
2974 | IF(A=='6') I=6 |
---|
2975 | IF(A=='7') I=7 |
---|
2976 | IF(A=='8') I=8 |
---|
2977 | IF(A=='9') I=9 |
---|
2978 | IF(A=='0') I=0 |
---|
2979 | if(i==-1) ndel=1 |
---|
2980 | IF(A=='a') I=1 |
---|
2981 | IF(A=='b') I=2 |
---|
2982 | IF(A=='c') I=3 |
---|
2983 | IF(A=='d') I=4 |
---|
2984 | IF(A=='e') I=5 |
---|
2985 | IF(A=='f') I=6 |
---|
2986 | IF(A=='g') I=7 |
---|
2987 | IF(A=='h') I=8 |
---|
2988 | IF(A=='i') I=9 |
---|
2989 | IF(A==' ') I=0 |
---|
2990 | IF(A=='o') I=0 |
---|
2991 | IF(A=='A') I=1 |
---|
2992 | IF(A=='B') I=2 |
---|
2993 | IF(A=='C') I=3 |
---|
2994 | IF(A=='D') I=4 |
---|
2995 | IF(A=='E') I=5 |
---|
2996 | IF(A=='F') I=6 |
---|
2997 | IF(A=='G') I=7 |
---|
2998 | IF(A=='H') I=8 |
---|
2999 | IF(A=='I') I=9 |
---|
3000 | IF(A=='O') I=0 |
---|
3001 | |
---|
3002 | |
---|
3003 | |
---|
3004 | END SUBROUTINE CHARINT |
---|
3005 | |
---|
3006 | |
---|
3007 | function check_j(j) |
---|
3008 | implicit none |
---|
3009 | integer check_j |
---|
3010 | INTEGER,INTENT(in),dimension(:)::j |
---|
3011 | integer i,no |
---|
3012 | IF(.NOT.C_%STABLE_DA) RETURN |
---|
3013 | |
---|
3014 | check_j=0 |
---|
3015 | |
---|
3016 | no=0 |
---|
3017 | do i=1,size(j) |
---|
3018 | no=j(i)+no |
---|
3019 | enddo |
---|
3020 | |
---|
3021 | if(no>c_%no) then |
---|
3022 | check_j=no |
---|
3023 | return |
---|
3024 | endif |
---|
3025 | |
---|
3026 | do i=c_%nv+1,size(j) |
---|
3027 | if(j(i)/=0) then |
---|
3028 | check_j=-i |
---|
3029 | endif |
---|
3030 | enddo |
---|
3031 | end function check_j |
---|
3032 | |
---|
3033 | |
---|
3034 | |
---|
3035 | function filter(j) |
---|
3036 | implicit none |
---|
3037 | real(dp) filter |
---|
3038 | integer i |
---|
3039 | integer,dimension(:)::j |
---|
3040 | |
---|
3041 | filter=1.0_dp |
---|
3042 | !do i=1,nd2+ndel |
---|
3043 | do i=1,nd2par |
---|
3044 | if(jfil(i)/=j(i)) filter=0.0_dp |
---|
3045 | enddo |
---|
3046 | |
---|
3047 | end function filter |
---|
3048 | |
---|
3049 | function filter_part(j) |
---|
3050 | implicit none |
---|
3051 | real(dp) filter_part |
---|
3052 | integer i |
---|
3053 | integer,dimension(:)::j |
---|
3054 | ! WRITE(6,*) jfilt(1:4) |
---|
3055 | ! WRITE(6,*)nd2part,nd2partt |
---|
3056 | filter_part=1.0_dp |
---|
3057 | !do i=1,nd2+ndel |
---|
3058 | do i=nd2part,nd2partt |
---|
3059 | if(jfilt(i)/=j(i)) filter_part=0.0_dp |
---|
3060 | enddo |
---|
3061 | |
---|
3062 | end function filter_part |
---|
3063 | |
---|
3064 | ! i/o routines |
---|
3065 | |
---|
3066 | SUBROUTINE pri(S1,MFILE,DEPS) |
---|
3067 | implicit none |
---|
3068 | INTEGER,INTENT(IN)::MFILE |
---|
3069 | REAL(DP),OPTIONAL,INTENT(INOUT)::DEPS |
---|
3070 | type (TAYLOR),INTENT(IN)::S1 |
---|
3071 | REAL(DP) PREC |
---|
3072 | |
---|
3073 | IF(PRESENT(DEPS)) THEN |
---|
3074 | PREC=-1.0_dp |
---|
3075 | CALL taylor_eps(PREC) |
---|
3076 | CALL taylor_eps(DEPS) |
---|
3077 | ENDIF |
---|
3078 | |
---|
3079 | ! if(old) then |
---|
3080 | if(print77) then |
---|
3081 | CALL DAPRI77(s1%i,MFILE) |
---|
3082 | else |
---|
3083 | CALL DAPRI(s1%i,MFILE) |
---|
3084 | endif |
---|
3085 | ! else |
---|
3086 | ! if(newprint) then |
---|
3087 | ! CALL newDAPRI(s1%j,MFILE) |
---|
3088 | ! else |
---|
3089 | ! if(print77) then |
---|
3090 | ! CALL oldDAPRI77(s1%j,MFILE) |
---|
3091 | ! else |
---|
3092 | ! CALL oldDAPRI(s1%j,MFILE) |
---|
3093 | ! endif |
---|
3094 | ! endif |
---|
3095 | ! endif |
---|
3096 | ! |
---|
3097 | IF(PRESENT(DEPS)) CALL taylor_eps(PREC) |
---|
3098 | |
---|
3099 | END SUBROUTINE pri |
---|
3100 | |
---|
3101 | SUBROUTINE REA(S1,MFILE) |
---|
3102 | implicit none |
---|
3103 | INTEGER,INTENT(in)::MFILE |
---|
3104 | type (TAYLOR),INTENT(IN)::S1 |
---|
3105 | |
---|
3106 | ! if(old) then |
---|
3107 | if(s1%i==0)call crap1("REA 1" ) ! call etall1(s1%i) |
---|
3108 | |
---|
3109 | if(read77) then |
---|
3110 | CALL DAREA77(s1%i,MFILE) |
---|
3111 | else |
---|
3112 | CALL DAREA(s1%i,MFILE) |
---|
3113 | endif |
---|
3114 | ! else |
---|
3115 | ! if(.NOT. ASSOCIATED(s1%j%r))call crap1("REA 2" ) ! call newetall(s1%j,1) |
---|
3116 | ! if(newread) then |
---|
3117 | ! CALL newDAREA(s1%j,MFILE) |
---|
3118 | ! else |
---|
3119 | ! if(read77) then |
---|
3120 | ! CALL oldDAREA77(s1%j,MFILE) |
---|
3121 | ! else |
---|
3122 | ! CALL oldDAREA(s1%j,MFILE) |
---|
3123 | ! endif |
---|
3124 | ! endif |
---|
3125 | ! endif |
---|
3126 | |
---|
3127 | END SUBROUTINE REA |
---|
3128 | |
---|
3129 | |
---|
3130 | ! Universal Taylor Routines (Sagan's Stuff) |
---|
3131 | |
---|
3132 | SUBROUTINE kill_uni(S2) |
---|
3133 | implicit none |
---|
3134 | type (UNIVERSAL_TAYLOR),INTENT(INOUT)::S2 |
---|
3135 | |
---|
3136 | DEALLOCATE(S2%N,S2%NV,S2%C,S2%J) |
---|
3137 | NULLIFY(S2%N,S2%NV,S2%C,S2%J) |
---|
3138 | |
---|
3139 | END SUBROUTINE kill_uni |
---|
3140 | |
---|
3141 | SUBROUTINE null_uni(S2,S1) |
---|
3142 | implicit none |
---|
3143 | type (UNIVERSAL_TAYLOR),INTENT(INOUT)::S2 |
---|
3144 | integer, intent(in):: s1 |
---|
3145 | IF(S1==0) THEN |
---|
3146 | NULLIFY(S2%N,S2%NV,S2%C,S2%J) |
---|
3147 | ELSEIF(S1==-1) THEN |
---|
3148 | DEALLOCATE(S2%N,S2%NV,S2%C,S2%J) |
---|
3149 | NULLIFY(S2%N,S2%NV,S2%C,S2%J) |
---|
3150 | ENDIF |
---|
3151 | END SUBROUTINE null_uni |
---|
3152 | |
---|
3153 | |
---|
3154 | SUBROUTINE ALLOC_U(S2,N,NV) |
---|
3155 | implicit none |
---|
3156 | type (UNIVERSAL_TAYLOR),INTENT(INOUT)::S2 |
---|
3157 | integer, intent(in):: N,NV |
---|
3158 | ALLOCATE(S2%N,S2%NV) |
---|
3159 | if(N==0) then |
---|
3160 | allocate(S2%C(1),S2%J(1,NV));S2%C(1)=0.0_dp;S2%J(:,:)=0; |
---|
3161 | else |
---|
3162 | allocate(S2%C(N),S2%J(N,NV)) |
---|
3163 | endif |
---|
3164 | S2%N=N |
---|
3165 | S2%NV=NV |
---|
3166 | END SUBROUTINE ALLOC_U |
---|
3167 | |
---|
3168 | SUBROUTINE fill_uni_r(S2,S1) !new sagan |
---|
3169 | implicit none |
---|
3170 | type (UNIVERSAL_TAYLOR),INTENT(INOUT)::S2 |
---|
3171 | real (dp), intent(in):: s1 |
---|
3172 | INTEGER n,J(LNV) |
---|
3173 | |
---|
3174 | |
---|
3175 | IF(ASSOCIATED(S2%N)) S2=-1 |
---|
3176 | S2=0 |
---|
3177 | CALL ALLOC_U(S2,1,nv) |
---|
3178 | J=0 |
---|
3179 | DO N=1,S2%NV |
---|
3180 | S2%J(1,N)=J(N) |
---|
3181 | ENDDO |
---|
3182 | S2%C(1)=S1 |
---|
3183 | |
---|
3184 | END SUBROUTINE fill_uni_r |
---|
3185 | |
---|
3186 | SUBROUTINE FILL_UNI(S2,S1) |
---|
3187 | implicit none |
---|
3188 | type (UNIVERSAL_TAYLOR),INTENT(INOUT)::S2 |
---|
3189 | type (TAYLOR), intent(in):: s1 |
---|
3190 | INTEGER ipresent,k,n,I,illa |
---|
3191 | real(dp) value |
---|
3192 | INTEGER, allocatable :: j(:) |
---|
3193 | call check_snake |
---|
3194 | |
---|
3195 | ! if(old) then |
---|
3196 | if(s1%i==0) call crap1("FILL_N 1") |
---|
3197 | ! else |
---|
3198 | ! IF (.NOT. ASSOCIATED(s1%j%r)) call crap1("FILL_N 2") |
---|
3199 | ! endif |
---|
3200 | |
---|
3201 | |
---|
3202 | IF(ASSOCIATED(S2%N)) S2=-1 |
---|
3203 | S2=0 |
---|
3204 | ipresent=1 |
---|
3205 | call dacycle(S1%I,ipresent,value,n) |
---|
3206 | CALL ALLOC_U(S2,N,c_%nv) |
---|
3207 | allocate(j(c_%nv)) |
---|
3208 | |
---|
3209 | do i=1,N |
---|
3210 | call dacycle(S1%I,i,value,illa,j) |
---|
3211 | S2%C(I)=value |
---|
3212 | DO k=1,S2%NV |
---|
3213 | S2%J(i,k)=J(k) |
---|
3214 | ENDDO |
---|
3215 | ENDDO |
---|
3216 | |
---|
3217 | deallocate(j) |
---|
3218 | |
---|
3219 | END SUBROUTINE FILL_UNI |
---|
3220 | |
---|
3221 | |
---|
3222 | |
---|
3223 | |
---|
3224 | SUBROUTINE REFILL_UNI(S1,S2) |
---|
3225 | implicit none |
---|
3226 | type (UNIVERSAL_TAYLOR),INTENT(IN)::S2 |
---|
3227 | type (TAYLOR), intent(inOUT):: s1 |
---|
3228 | INTEGER I,K,J(LNV) |
---|
3229 | logical(lp) DOIT |
---|
3230 | |
---|
3231 | ! if(old) then |
---|
3232 | if(s1%i==0) call crap1("REFILL_N 1") |
---|
3233 | ! else |
---|
3234 | ! IF (.NOT. ASSOCIATED(s1%j%r)) call crap1("REFILL_N 2") |
---|
3235 | ! endif |
---|
3236 | |
---|
3237 | |
---|
3238 | S1=0.0_dp |
---|
3239 | |
---|
3240 | IF(.not.ASSOCIATED(S2%N)) THEN |
---|
3241 | w_p=0 |
---|
3242 | w_p%nc=1 |
---|
3243 | w_p%fc='(1((1X,A72),/))' |
---|
3244 | w_p%c(1)=" ERROR IN REFILL_N: UNIVERSAL_TAYLOR DOES NOT EXIST" |
---|
3245 | ! call !write_e(123) |
---|
3246 | ENDIF |
---|
3247 | J=0 |
---|
3248 | DO I=1,S2%N |
---|
3249 | DOIT=.TRUE. |
---|
3250 | IF(S2%NV>NV) THEN |
---|
3251 | K=NV |
---|
3252 | DO WHILE(DOIT.AND.K<=S2%NV) |
---|
3253 | IF(S2%J(I,K)/=0) DOIT=.FALSE. |
---|
3254 | K=K+1 |
---|
3255 | ENDDO |
---|
3256 | ENDIF |
---|
3257 | |
---|
3258 | IF(DOIT) THEN |
---|
3259 | DO K=1,NV |
---|
3260 | J(K)=S2%J(I,K) |
---|
3261 | ENDDO |
---|
3262 | CALL POK(S1,J,S2%C(I)) |
---|
3263 | ENDIF |
---|
3264 | ENDDO |
---|
3265 | |
---|
3266 | END SUBROUTINE REFILL_UNI |
---|
3267 | |
---|
3268 | |
---|
3269 | !_________________________________________________________________________________ |
---|
3270 | |
---|
3271 | |
---|
3272 | subroutine printunitaylor(ut,iunit) |
---|
3273 | implicit none |
---|
3274 | type(universal_taylor) :: ut |
---|
3275 | integer :: iunit |
---|
3276 | integer :: i,ii |
---|
3277 | |
---|
3278 | if (.not. associated(ut%n)) then |
---|
3279 | write(iunit,'(A)') ' UNIVERSAL_TAYLOR IS EMPTY (NOT ASSOCIATED)' |
---|
3280 | write(6,'(A)') ' UNIVERSAL_TAYLOR IS EMPTY (NOT ASSOCIATED)' |
---|
3281 | return |
---|
3282 | endif |
---|
3283 | |
---|
3284 | write(iunit,'(/1X,A,I5,A,I5,A/1X,A/)') 'UNIV_TAYLOR NO =',ut%n,', NV =',ut%nv,', INA = unita',& |
---|
3285 | '*********************************************' |
---|
3286 | if(ut%n /= 0) then |
---|
3287 | write(iunit,'(A)') ' I COEFFICIENT ORDER EXPONENTS' |
---|
3288 | else |
---|
3289 | write(iunit,'(A)') ' ALL COMPONENTS 0.0_dp ' |
---|
3290 | endif |
---|
3291 | |
---|
3292 | do i = 1,ut%n |
---|
3293 | write(iunit,'(I6,2X,G21.14,I5,4X,18(2I2,1X))') i,ut%c(i),sum(ut%j(i,:)),(ut%j(i,ii),ii=1,ut%nv) |
---|
3294 | if( .not. print77) then |
---|
3295 | write(iunit,*) ut%c(i) |
---|
3296 | endif |
---|
3297 | enddo |
---|
3298 | |
---|
3299 | write(iunit,'(A)') ' ' |
---|
3300 | |
---|
3301 | end subroutine printunitaylor |
---|
3302 | |
---|
3303 | |
---|
3304 | ! End of Universal Taylor Routines |
---|
3305 | |
---|
3306 | |
---|
3307 | |
---|
3308 | |
---|
3309 | ! Warning Routines |
---|
3310 | |
---|
3311 | subroutine crap1(STRING) |
---|
3312 | implicit none |
---|
3313 | CHARACTER(*) STRING |
---|
3314 | |
---|
3315 | w_p=0 |
---|
3316 | w_p%nc=2 |
---|
3317 | w_p%fc='((1X,A72,/),(1X,A72))' |
---|
3318 | w_p%c(1)= "ERROR IN :" |
---|
3319 | w_p%c(2)= STRING |
---|
3320 | ! call !write_e(3478) |
---|
3321 | |
---|
3322 | end subroutine crap1 |
---|
3323 | |
---|
3324 | SUBROUTINE real_stop() |
---|
3325 | implicit none |
---|
3326 | integer i(1),j |
---|
3327 | |
---|
3328 | w_p=0 |
---|
3329 | w_p%nc=3 |
---|
3330 | write(6,*) " You are using a kind(1.0_dp) " |
---|
3331 | write(6,*)" set real_warning to false to permit this " |
---|
3332 | write(6,*)" write 1 to continue or -1 for a crash " |
---|
3333 | call read(j) |
---|
3334 | i(j)=0 |
---|
3335 | real_warning=.false. |
---|
3336 | |
---|
3337 | END SUBROUTINE real_stop |
---|
3338 | |
---|
3339 | |
---|
3340 | SUBROUTINE ndum_warning_user() |
---|
3341 | implicit none |
---|
3342 | integer ipause,II(0:1) |
---|
3343 | |
---|
3344 | |
---|
3345 | w_p=0 |
---|
3346 | w_p%nc=3 |
---|
3347 | w_p%fc='(3((1X,A72),/))' |
---|
3348 | w_p%c(1)= " *****************************************************************" |
---|
3349 | w_p%c(2)= " * Should never be here in New Linked List Scheme *" |
---|
3350 | w_p%c(3)= " *****************************************************************" |
---|
3351 | w_p=0 |
---|
3352 | w_p%nc=1 |
---|
3353 | w_p%fc='(1(1X,A72),/))' |
---|
3354 | w_p%c(1)= " do you want a crash? " |
---|
3355 | ! call !write_e |
---|
3356 | call read(ipause) |
---|
3357 | ii(2000*ipause)=0 |
---|
3358 | |
---|
3359 | end SUBROUTINE ndum_warning_user |
---|
3360 | |
---|
3361 | ! End of Warning Routines |
---|
3362 | |
---|
3363 | ! linked list of da for scratch levels |
---|
3364 | |
---|
3365 | SUBROUTINE Set_Up( L ) ! Sets up a layout: gives a unique negative index |
---|
3366 | implicit none |
---|
3367 | TYPE (dalevel) L |
---|
3368 | call null_it(L) |
---|
3369 | ALLOCATE(L%n); |
---|
3370 | ALLOCATE(L%CLOSED); |
---|
3371 | L%closed=.FALSE. |
---|
3372 | L%N=0 |
---|
3373 | END SUBROUTINE Set_Up |
---|
3374 | |
---|
3375 | SUBROUTINE de_Set_Up( L ) ! deallocates layout content |
---|
3376 | implicit none |
---|
3377 | TYPE (dalevel) L |
---|
3378 | deallocate(L%closed); |
---|
3379 | deallocate(L%n); |
---|
3380 | END SUBROUTINE de_Set_Up |
---|
3381 | |
---|
3382 | |
---|
3383 | |
---|
3384 | SUBROUTINE null_it( L ) ! Nullifies layout content |
---|
3385 | implicit none |
---|
3386 | TYPE (dalevel), intent(inout) :: L |
---|
3387 | nullify(L%N ) |
---|
3388 | nullify(L%CLOSED ) |
---|
3389 | nullify(L%PRESENT ) |
---|
3390 | ! |
---|
3391 | nullify(L%END ) |
---|
3392 | nullify(L%START ) |
---|
3393 | nullify(L%START_GROUND )! STORE THE GROUNDED VALUE OF START DURING CIRCULAR SCANNING |
---|
3394 | nullify(L%END_GROUND )! STORE THE GROUNDED VALUE OF END DURING CIRCULAR SCANNING |
---|
3395 | END SUBROUTINE null_it |
---|
3396 | |
---|
3397 | SUBROUTINE LINE_L(L,doneit) ! makes into line temporarily |
---|
3398 | implicit none |
---|
3399 | TYPE (DALEVEL) L |
---|
3400 | logical(lp) doneit |
---|
3401 | doneit=.false. |
---|
3402 | if(L%closed) then |
---|
3403 | if(associated(L%end%next)) then |
---|
3404 | L%end%next=>L%start_ground |
---|
3405 | doneit=.true. |
---|
3406 | endif |
---|
3407 | if(associated(L%start%previous)) then |
---|
3408 | L%start%previous=>L%end_ground |
---|
3409 | endif |
---|
3410 | endif |
---|
3411 | END SUBROUTINE LINE_L |
---|
3412 | |
---|
3413 | SUBROUTINE RING_L(L,doit) ! Brings back to ring if needed |
---|
3414 | implicit none |
---|
3415 | TYPE (DALEVEL) L |
---|
3416 | logical(lp) doit |
---|
3417 | if(L%closed.and.doit) then |
---|
3418 | if(.NOT.(associated(L%end%next))) then |
---|
3419 | L%start_ground=>L%end%next ! saving grounded pointer |
---|
3420 | L%end%next=>L%start |
---|
3421 | endif |
---|
3422 | if(.NOT.(associated(L%start%previous))) then |
---|
3423 | L%end_ground=>L%start%previous ! saving grounded pointer |
---|
3424 | L%start%previous=>L%end |
---|
3425 | endif |
---|
3426 | endif |
---|
3427 | END SUBROUTINE RING_L |
---|
3428 | |
---|
3429 | SUBROUTINE APPEND_DA( L ) ! Standard append that clones everything |
---|
3430 | implicit none |
---|
3431 | TYPE (dascratch), POINTER :: Current |
---|
3432 | TYPE (DALEVEL), TARGET,intent(inout):: L |
---|
3433 | logical(lp) doneit |
---|
3434 | CALL LINE_L(L,doneit) |
---|
3435 | L%N=L%N+1 |
---|
3436 | nullify(current);ALLOCATE(current); |
---|
3437 | |
---|
3438 | call alloc_DA(current) |
---|
3439 | |
---|
3440 | if(L%N==1) current%next=> L%start |
---|
3441 | Current % previous => L % end ! point it to next fibre |
---|
3442 | if(L%N>1) THEN |
---|
3443 | L % end % next => current ! |
---|
3444 | ENDIF |
---|
3445 | |
---|
3446 | L % end => Current |
---|
3447 | if(L%N==1) L%start=> Current |
---|
3448 | L%PRESENT=>CURRENT ! ALWAYS IF APPENDING |
---|
3449 | CALL RING_L(L,doneit) |
---|
3450 | END SUBROUTINE APPEND_DA |
---|
3451 | |
---|
3452 | SUBROUTINE INSERT_DA( L ) ! Standard append that clones everything |
---|
3453 | implicit none |
---|
3454 | logical(lp) :: doneitt=.true. |
---|
3455 | TYPE (dascratch), POINTER :: Current |
---|
3456 | TYPE (DALEVEL), TARGET,intent(inout):: L |
---|
3457 | IF(L%N>1.AND.(.NOT.ASSOCIATED(L%PRESENT,L%END))) THEN |
---|
3458 | |
---|
3459 | L%N=L%N+1 |
---|
3460 | nullify(current);ALLOCATE(current); |
---|
3461 | |
---|
3462 | call alloc_DA(current) |
---|
3463 | |
---|
3464 | Current % previous => L % PRESENT ! 2P -> 2 |
---|
3465 | Current % NEXT => L % PRESENT%NEXT ! 2P -> 3 |
---|
3466 | L%PRESENT%NEXT=> CURRENT ! 2 -> 2P |
---|
3467 | Current % NEXT%PREVIOUS => CURRENT ! 3 -> 2P |
---|
3468 | L%PRESENT=>CURRENT ! 2P BECOMES 3 |
---|
3469 | |
---|
3470 | ELSE |
---|
3471 | |
---|
3472 | CALL APPEND_DA( L ) |
---|
3473 | if(L%N==1) THEN |
---|
3474 | L%CLOSED=.TRUE. |
---|
3475 | CALL RING_L(L,doneitt) |
---|
3476 | ENDIF |
---|
3477 | |
---|
3478 | ENDIF |
---|
3479 | END SUBROUTINE INSERT_DA |
---|
3480 | |
---|
3481 | SUBROUTINE alloc_DA( c ) ! Does the full allocation of fibre and initialization of internal variables |
---|
3482 | implicit none |
---|
3483 | type(dascratch),pointer:: c |
---|
3484 | ALLOCATE(C%T) |
---|
3485 | CALL ALLOC(C%T) |
---|
3486 | NULLIFY(C%NEXT) |
---|
3487 | NULLIFY(C%PREVIOUS) |
---|
3488 | |
---|
3489 | end SUBROUTINE alloc_DA |
---|
3490 | |
---|
3491 | SUBROUTINE kill_DALEVEL( L ) ! Destroys a layout |
---|
3492 | implicit none |
---|
3493 | TYPE (DASCRATCH), POINTER :: Current |
---|
3494 | TYPE (DALEVEL) L |
---|
3495 | logical(lp) doneit |
---|
3496 | CALL LINE_L(L,doneit) |
---|
3497 | nullify(current) |
---|
3498 | Current => L % end ! end at the end |
---|
3499 | DO WHILE (ASSOCIATED(L % end)) |
---|
3500 | L % end => Current % previous ! update the end before disposing |
---|
3501 | call dealloc_DASCRATCH(Current) |
---|
3502 | Current => L % end ! alias of last fibre again |
---|
3503 | L%N=L%N-1 |
---|
3504 | END DO |
---|
3505 | call de_set_up(L) |
---|
3506 | END SUBROUTINE kill_DALEVEL |
---|
3507 | |
---|
3508 | SUBROUTINE dealloc_DASCRATCH( c ) ! destroys internal data if it is not pointing (i.e. not a parent) |
---|
3509 | implicit none |
---|
3510 | type(DASCRATCH),pointer :: c |
---|
3511 | IF(ASSOCIATED(C)) THEN |
---|
3512 | CALL KILL(C%T) |
---|
3513 | IF(ASSOCIATED(C%T)) DEALLOCATE(C%T) |
---|
3514 | ! IF(ASSOCIATED(C%NEXT)) DEALLOCATE(C%NEXT) |
---|
3515 | ! IF(ASSOCIATED(C%PREVIOUS)) DEALLOCATE(C%PREVIOUS) |
---|
3516 | deallocate(c); |
---|
3517 | ENDIF |
---|
3518 | end SUBROUTINE dealloc_DASCRATCH |
---|
3519 | |
---|
3520 | SUBROUTINE set_up_level() |
---|
3521 | implicit none |
---|
3522 | integer i |
---|
3523 | do i=1,ndumt |
---|
3524 | call set_up(scratchda(i)) |
---|
3525 | ! do j=1,n |
---|
3526 | ! call INSERT_da(scratchda(i)) |
---|
3527 | ! enddo |
---|
3528 | ! scratchda(i)%CLOSED=.TRUE. |
---|
3529 | ! CALL RING_L(scratchda(i),.TRUE.) |
---|
3530 | enddo |
---|
3531 | |
---|
3532 | end SUBROUTINE set_up_level |
---|
3533 | |
---|
3534 | SUBROUTINE report_level() |
---|
3535 | implicit none |
---|
3536 | integer i |
---|
3537 | if(associated(scratchda(1)%n)) then |
---|
3538 | do i=1,ndumt |
---|
3539 | w_p=0 |
---|
3540 | w_p%nc=1 |
---|
3541 | w_p%fc='(1((1X,A72)))' |
---|
3542 | write(6,'(a6,1x,i4,a5,1x,i4,1x,a7)') "Level ",i, " has ",scratchda(i)%n, "Taylors" |
---|
3543 | ! write(w_p%c(1),'(a6,1x,i4,a5,1x,i4,1x,a7)') "Level ",i, " has ",scratchda(i)%n, "Taylors" |
---|
3544 | ! ! call !write_e |
---|
3545 | enddo |
---|
3546 | endif |
---|
3547 | END SUBROUTINE report_level |
---|
3548 | |
---|
3549 | ! end linked list of da for scratch levels |
---|
3550 | |
---|
3551 | ! Assignments Routines |
---|
3552 | |
---|
3553 | subroutine ASSIGN() |
---|
3554 | implicit none |
---|
3555 | integer i |
---|
3556 | do i=1,ndumt |
---|
3557 | iassdoluser(i)=0 |
---|
3558 | iass0user(i)=0 |
---|
3559 | enddo |
---|
3560 | ! if(old) then |
---|
3561 | CALL ETALL1(DUMMY) |
---|
3562 | call etall1(temp) |
---|
3563 | ! else |
---|
3564 | ! CALL allocnewda(DUMMYl) |
---|
3565 | ! call allocnewda(templ) |
---|
3566 | ! endif |
---|
3567 | CALL set_up_level |
---|
3568 | end subroutine ASSIGN |
---|
3569 | |
---|
3570 | subroutine DEASSIGN() |
---|
3571 | implicit none |
---|
3572 | integer i |
---|
3573 | do i=1,ndumt |
---|
3574 | iassdoluser(i)=0 |
---|
3575 | iass0user(i)=0 |
---|
3576 | enddo |
---|
3577 | ! if(old) then |
---|
3578 | CALL DADAL1(DUMMY) |
---|
3579 | call DADAL1(temp) |
---|
3580 | ! else |
---|
3581 | ! CALL KILLnewdaS(DUMMYl) |
---|
3582 | ! call KILLnewdaS(templ) |
---|
3583 | ! endif |
---|
3584 | do i=1,ndumt |
---|
3585 | CALL kill_DALEVEL(scratchda(I)) |
---|
3586 | ENDDO |
---|
3587 | end subroutine DEASSIGN |
---|
3588 | |
---|
3589 | subroutine ASStaylor(s1) |
---|
3590 | implicit none |
---|
3591 | TYPE (taylor) s1 |
---|
3592 | ! lastmaster=master ! 2002.12.13 |
---|
3593 | |
---|
3594 | select case(master) |
---|
3595 | case(0:ndumt-1) |
---|
3596 | master=master+1 |
---|
3597 | case(ndumt) |
---|
3598 | write(6,*) " cannot indent anymore ",ndumt |
---|
3599 | w_p=0 |
---|
3600 | w_p%nc=1 |
---|
3601 | w_p=(/" cannot indent anymore "/) |
---|
3602 | w_p%fc='(1((1X,A72),/))' |
---|
3603 | ! call !write_e(100) |
---|
3604 | master=sqrt(-dble(master)) |
---|
3605 | end select |
---|
3606 | ! write(26,*) " taylor ",master |
---|
3607 | call ass0(s1) |
---|
3608 | |
---|
3609 | end subroutine ASStaylor |
---|
3610 | |
---|
3611 | subroutine ass0(s1) |
---|
3612 | implicit none |
---|
3613 | integer ipause, mypause |
---|
3614 | TYPE (taylor) s1 |
---|
3615 | |
---|
3616 | IF(MASTER>NDUMT.or.master==0) THEN |
---|
3617 | WRITE(6,*) "more scratch level needed ",master,NDUMT |
---|
3618 | ipause=mypause(123) |
---|
3619 | write(6,*) 1/sqrt(-dble(1000+master)) |
---|
3620 | stop 123 |
---|
3621 | ENDIF |
---|
3622 | |
---|
3623 | if(.not.no_ndum_check) iass0user(master)=iass0user(master)+1 |
---|
3624 | if(iass0user(master)>scratchda(master)%n) then |
---|
3625 | call INSERT_DA( scratchda(master) ) |
---|
3626 | ELSE |
---|
3627 | scratchda(master)%PRESENT=>scratchda(master)%PRESENT%NEXT |
---|
3628 | ENDIF |
---|
3629 | ! if(old) then |
---|
3630 | s1%i=scratchda(master)%PRESENT%T%i |
---|
3631 | ! else |
---|
3632 | ! s1%j=scratchda(master)%PRESENT%T%j |
---|
3633 | ! endif |
---|
3634 | |
---|
3635 | |
---|
3636 | end subroutine ASS0 |
---|
3637 | |
---|
3638 | |
---|
3639 | ! remove small numbers |
---|
3640 | |
---|
3641 | SUBROUTINE clean_taylor(S1,S2,prec) |
---|
3642 | implicit none |
---|
3643 | type (TAYLOR),INTENT(INOUT)::S2 |
---|
3644 | type (TAYLOR), intent(INOUT):: s1 |
---|
3645 | real(dp) prec |
---|
3646 | INTEGER ipresent,k,n,I,illa |
---|
3647 | real(dp) value |
---|
3648 | INTEGER, allocatable :: j(:) |
---|
3649 | type (TAYLOR) t |
---|
3650 | |
---|
3651 | call alloc(t) |
---|
3652 | t=0.0_dp |
---|
3653 | ipresent=1 |
---|
3654 | call dacycle(S1%I,ipresent,value,n) |
---|
3655 | |
---|
3656 | allocate(j(c_%nv)) |
---|
3657 | |
---|
3658 | do i=1,N |
---|
3659 | call dacycle(S1%I,i,value,illa,j) |
---|
3660 | if(abs(value)>prec) then |
---|
3661 | t=t+(value.mono.j) |
---|
3662 | endif |
---|
3663 | ENDDO |
---|
3664 | s2=t |
---|
3665 | deallocate(j) |
---|
3666 | call kill(t) |
---|
3667 | |
---|
3668 | END SUBROUTINE clean_taylor |
---|
3669 | |
---|
3670 | |
---|
3671 | SUBROUTINE clean_pbfield(S1,S2,prec) |
---|
3672 | implicit none |
---|
3673 | type (pbfield),INTENT(INOUT)::S2 |
---|
3674 | type (pbfield), intent(INOUT):: s1 |
---|
3675 | real(dp) prec |
---|
3676 | |
---|
3677 | call clean_taylor(s1%h,s2%h,prec) |
---|
3678 | |
---|
3679 | END SUBROUTINE clean_pbfield |
---|
3680 | |
---|
3681 | SUBROUTINE clean_pbresonance (S1,S2,prec) |
---|
3682 | implicit none |
---|
3683 | type (pbresonance),INTENT(INOUT)::S2 |
---|
3684 | type (pbresonance), intent(INOUT):: s1 |
---|
3685 | real(dp) prec |
---|
3686 | |
---|
3687 | call clean_pbfield(s1%cos,s2%cos,prec) |
---|
3688 | call clean_pbfield(s1%sin,s2%sin,prec) |
---|
3689 | |
---|
3690 | END SUBROUTINE clean_pbresonance |
---|
3691 | |
---|
3692 | SUBROUTINE clean_damap(S1,S2,prec) |
---|
3693 | implicit none |
---|
3694 | type (damap),INTENT(INOUT)::S2 |
---|
3695 | type (damap), intent(INOUT):: s1 |
---|
3696 | real(dp) prec |
---|
3697 | integer i |
---|
3698 | |
---|
3699 | do i=1,c_%nd2 |
---|
3700 | call clean_taylor(s1%v(i),s2%v(i),prec) |
---|
3701 | enddo |
---|
3702 | |
---|
3703 | |
---|
3704 | END SUBROUTINE clean_damap |
---|
3705 | |
---|
3706 | SUBROUTINE clean_vecfield(S1,S2,prec) |
---|
3707 | implicit none |
---|
3708 | type (vecfield),INTENT(INOUT)::S2 |
---|
3709 | type (vecfield), intent(INOUT):: s1 |
---|
3710 | real(dp) prec |
---|
3711 | integer i |
---|
3712 | |
---|
3713 | do i=1,c_%nd2 |
---|
3714 | call clean_taylor(s1%v(i),s2%v(i),prec) |
---|
3715 | enddo |
---|
3716 | |
---|
3717 | |
---|
3718 | END SUBROUTINE clean_vecfield |
---|
3719 | |
---|
3720 | SUBROUTINE clean_vecresonance(S1,S2,prec) |
---|
3721 | implicit none |
---|
3722 | type (vecresonance),INTENT(INOUT)::S2 |
---|
3723 | type (vecresonance), intent(INOUT):: s1 |
---|
3724 | real(dp) prec |
---|
3725 | integer i |
---|
3726 | |
---|
3727 | |
---|
3728 | call clean_vecfield(s1%cos,s2%cos,prec) |
---|
3729 | call clean_vecfield(s1%sin,s2%sin,prec) |
---|
3730 | |
---|
3731 | |
---|
3732 | |
---|
3733 | END SUBROUTINE clean_vecresonance |
---|
3734 | |
---|
3735 | SUBROUTINE clean_onelieexponent(S1,S2,prec) |
---|
3736 | implicit none |
---|
3737 | type (onelieexponent),INTENT(INOUT)::S2 |
---|
3738 | type (onelieexponent), intent(INOUT):: s1 |
---|
3739 | real(dp) prec |
---|
3740 | integer i |
---|
3741 | |
---|
3742 | |
---|
3743 | call clean_vecfield(s1%vector,s2%vector,prec) |
---|
3744 | call clean_pbfield(s1%pb,s2%pb,prec) |
---|
3745 | |
---|
3746 | |
---|
3747 | |
---|
3748 | END SUBROUTINE clean_onelieexponent |
---|
3749 | |
---|
3750 | |
---|
3751 | ! remove small numbers |
---|
3752 | |
---|
3753 | SUBROUTINE clean_complextaylor(S1,S2,prec) |
---|
3754 | implicit none |
---|
3755 | type (complextaylor),INTENT(INOUT)::S2 |
---|
3756 | type (complextaylor), intent(INOUT):: s1 |
---|
3757 | real(dp) prec |
---|
3758 | |
---|
3759 | call clean_taylor(S1%r,S1%r,prec) |
---|
3760 | call clean_taylor(S1%i,S1%i,prec) |
---|
3761 | |
---|
3762 | |
---|
3763 | END SUBROUTINE clean_complextaylor |
---|
3764 | |
---|
3765 | SUBROUTINE clean_gmap(S1,s2,prec) |
---|
3766 | implicit none |
---|
3767 | type (gmap),INTENT(INOUT)::S1 |
---|
3768 | type (gmap),INTENT(INOUT)::S2 |
---|
3769 | real(dp) prec |
---|
3770 | INTEGER I |
---|
3771 | |
---|
3772 | DO I=1,s1%n |
---|
3773 | CALL clean_taylor(S1%V(I),S2%V(I),prec) |
---|
3774 | ENDDO |
---|
3775 | |
---|
3776 | END SUBROUTINE clean_gmap |
---|
3777 | |
---|
3778 | |
---|
3779 | |
---|
3780 | |
---|
3781 | END MODULE tpsa |
---|